File: clsyvt.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 (967 lines) | stat: -rw-r--r-- 31,514 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
!-------------------------------------------------------------------------------

! 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 clsyvt.f90
!>
!> \brief Symmetry boundary conditions for vectors and tensors.
!>
!> Correspond to the code icodcl(ivar) = 4.
!-------------------------------------------------------------------------------

!-------------------------------------------------------------------------------
! Arguments
!______________________________________________________________________________.
!  mode           name          role                                           !
!______________________________________________________________________________!
!> \param[in]     nscal         total number of scalars
!> \param[in,out] icodcl        face boundary condition code:
!>                               - 1 Dirichlet
!>                               - 3 Neumann
!>                               - 4 sliding and
!>                                 \f$ \vect{u} \cdot \vect{n} = 0 \f$
!>                               - 5 smooth wall and
!>                                 \f$ \vect{u} \cdot \vect{n} = 0 \f$
!>                               - 6 rough wall and
!>                                 \f$ \vect{u} \cdot \vect{n} = 0 \f$
!>                               - 9 free inlet/outlet
!>                                 (input mass flux blocked to 0)
!> \param[in,out] rcodcl        boundary condition values:
!>                               - rcodcl(1) value of the dirichlet
!>                               - rcodcl(2) value of the exterior exchange
!>                                 coefficient (infinite if no exchange)
!>                               - rcodcl(3) value flux density
!>                                 (negative if gain) in w/m2 or roughness
!>                                 in m if icodcl=6
!>                                 -# for the velocity \f$ (\mu+\mu_T)
!>                                    \gradv \vect{u} \cdot \vect{n}  \f$
!>                                 -# for the pressure \f$ \Delta t
!>                                    \grad P \cdot \vect{n}  \f$
!>                                 -# for a scalar \f$ cp \left( K +
!>                                     \dfrac{K_T}{\sigma_T} \right)
!>                                     \grad T \cdot \vect{n} \f$
!> \param[in]     velipb        value of the velocity at \f$ \centip \f$
!>                               of boundary cells
!> \param[in]     rijipb        value of \f$ R_{ij} \f$ at \f$ \centip \f$
!>                               of boundary cells
!_______________________________________________________________________________


subroutine clsyvt &
 ( nscal  , icodcl ,                                     &
   rcodcl ,                                              &
   velipb , rijipb )

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

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

use paramx
use numvar
use optcal
use cstphy
use cstnum
use pointe
use entsor
use albase
use field
use mesh

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

implicit none

! Arguments

integer          nscal

integer          icodcl(nfabor,nvarcl)

double precision rcodcl(nfabor,nvarcl,3)
double precision velipb(nfabor,ndim), rijipb(nfabor,6)

! Local variables

integer          ifac, ii, isou
integer          iel
integer          iscal , ivar

double precision rnx, rny, rnz, rxnn
double precision upx, upy, upz, usn
double precision tx, ty, tz, txn, t2x, t2y, t2z
double precision clsyme
double precision eloglo(3,3), alpha(6,6)
double precision srfbnf, rcodcn, hint, visclc, visctc, distbf
double precision cpp
double precision vismsh(3), hintv(3)
double precision visci(3,3), fikis, viscis, distfi
double precision fcoefa(6), fcoefb(6), fcofaf(6), fcofbf(6), fcofad(6), fcofbd(6)

double precision, dimension(:,:), pointer :: coefau, cofafu, cfaale, claale
double precision, dimension(:,:), pointer :: visten
double precision, dimension(:,:,:), pointer :: coefbu, cofbfu, cfbale, clbale
double precision, dimension(:), pointer :: coefa_r11, coefaf_r11, coefad_r11
double precision, dimension(:), pointer :: coefb_r11, coefbf_r11, coefbd_r11
double precision, dimension(:), pointer :: coefa_r22, coefaf_r22, coefad_r22
double precision, dimension(:), pointer :: coefb_r22, coefbf_r22, coefbd_r22
double precision, dimension(:), pointer :: coefa_r33, coefaf_r33, coefad_r33
double precision, dimension(:), pointer :: coefb_r33, coefbf_r33, coefbd_r33
double precision, dimension(:), pointer :: coefa_r12, coefaf_r12, coefad_r12
double precision, dimension(:), pointer :: coefb_r12, coefbf_r12, coefbd_r12
double precision, dimension(:), pointer :: coefa_r13, coefaf_r13, coefad_r13
double precision, dimension(:), pointer :: coefb_r13, coefbf_r13, coefbd_r13
double precision, dimension(:), pointer :: coefa_r23, coefaf_r23, coefad_r23
double precision, dimension(:), pointer :: coefb_r23, coefbf_r23, coefbd_r23
double precision, dimension(:,:), pointer :: coefa_rij, coefaf_rij, coefad_rij
double precision, dimension(:,:,:), pointer :: coefb_rij, coefbf_rij, coefbd_rij

double precision, dimension(:), pointer :: viscl, visct
double precision, dimension(:), pointer :: cpro_visma1, cpro_visma2, cpro_visma3

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

!===============================================================================
! 1. Initializations
!===============================================================================

cpp = 0.d0

if (itytur.eq.3 .and. idirsm.eq.1) call field_get_val_v(ivsten, visten)

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

! Boundary Conditions

call field_get_coefa_v(ivarfl(iu), coefau)
call field_get_coefb_v(ivarfl(iu), coefbu)
call field_get_coefaf_v(ivarfl(iu), cofafu)
call field_get_coefbf_v(ivarfl(iu), cofbfu)

if (itytur.eq.3) then
  if (irijco.eq.1) then
    call field_get_coefa_v(ivarfl(irij), coefa_rij)
    call field_get_coefb_v(ivarfl(irij), coefb_rij)
    call field_get_coefaf_v(ivarfl(irij), coefaf_rij)
    call field_get_coefbf_v(ivarfl(irij), coefbf_rij)
    call field_get_coefad_v(ivarfl(irij), coefad_rij)
    call field_get_coefbd_v(ivarfl(irij), coefbd_rij)

    coefb_r11 => null()
    coefaf_r11 => null()
    coefbf_r11 => null()
    coefad_r11 => null()
    coefbd_r11 => null()

    coefa_r22 => null()
    coefb_r22 => null()
    coefaf_r22 => null()
    coefbf_r22 => null()
    coefad_r22 => null()
    coefbd_r22 => null()

    coefa_r33 => null()
    coefb_r33 => null()
    coefaf_r33 => null()
    coefbf_r33 => null()
    coefad_r33 => null()
    coefbd_r33 => null()

    coefa_r12 => null()
    coefb_r12 => null()
    coefaf_r12 => null()
    coefbf_r12 => null()
    coefad_r12 => null()
    coefbd_r12 => null()

    coefa_r13 => null()
    coefb_r13 => null()
    coefaf_r13 => null()
    coefbf_r13 => null()
    coefad_r13 => null()
    coefbd_r13 => null()

    coefa_r23 => null()
    coefb_r23 => null()
    coefaf_r23 => null()
    coefbf_r23 => null()
    coefad_r23 => null()
    coefbd_r23 => null()


  else

    call field_get_coefa_s(ivarfl(ir11), coefa_r11)
    call field_get_coefb_s(ivarfl(ir11), coefb_r11)
    call field_get_coefaf_s(ivarfl(ir11), coefaf_r11)
    call field_get_coefbf_s(ivarfl(ir11), coefbf_r11)
    call field_get_coefad_s(ivarfl(ir11), coefad_r11)
    call field_get_coefbd_s(ivarfl(ir11), coefbd_r11)

    call field_get_coefa_s(ivarfl(ir22), coefa_r22)
    call field_get_coefb_s(ivarfl(ir22), coefb_r22)
    call field_get_coefaf_s(ivarfl(ir22), coefaf_r22)
    call field_get_coefbf_s(ivarfl(ir22), coefbf_r22)
    call field_get_coefad_s(ivarfl(ir22), coefad_r22)
    call field_get_coefbd_s(ivarfl(ir22), coefbd_r22)

    call field_get_coefa_s(ivarfl(ir33), coefa_r33)
    call field_get_coefb_s(ivarfl(ir33), coefb_r33)
    call field_get_coefaf_s(ivarfl(ir33), coefaf_r33)
    call field_get_coefbf_s(ivarfl(ir33), coefbf_r33)
    call field_get_coefad_s(ivarfl(ir33), coefad_r33)
    call field_get_coefbd_s(ivarfl(ir33), coefbd_r33)

    call field_get_coefa_s(ivarfl(ir12), coefa_r12)
    call field_get_coefb_s(ivarfl(ir12), coefb_r12)
    call field_get_coefaf_s(ivarfl(ir12), coefaf_r12)
    call field_get_coefbf_s(ivarfl(ir12), coefbf_r12)
    call field_get_coefad_s(ivarfl(ir12), coefad_r12)
    call field_get_coefbd_s(ivarfl(ir12), coefbd_r12)

    call field_get_coefa_s(ivarfl(ir13), coefa_r13)
    call field_get_coefb_s(ivarfl(ir13), coefb_r13)
    call field_get_coefaf_s(ivarfl(ir13), coefaf_r13)
    call field_get_coefbf_s(ivarfl(ir13), coefbf_r13)
    call field_get_coefad_s(ivarfl(ir13), coefad_r13)
    call field_get_coefbd_s(ivarfl(ir13), coefbd_r13)

    call field_get_coefa_s(ivarfl(ir23), coefa_r23)
    call field_get_coefb_s(ivarfl(ir23), coefb_r23)
    call field_get_coefaf_s(ivarfl(ir23), coefaf_r23)
    call field_get_coefbf_s(ivarfl(ir23), coefbf_r23)
    call field_get_coefad_s(ivarfl(ir23), coefad_r23)
    call field_get_coefbd_s(ivarfl(ir23), coefbd_r23)

    coefa_rij => null()
    coefb_rij => null()
    coefaf_rij => null()
    coefbf_rij => null()
    coefad_rij => null()
    coefbd_rij => null()
  endif
endif

! --- Begin the loop over boundary faces
do ifac = 1, nfabor

  ! --- Test sur la presence d'une condition de symetrie vitesse : debut
  if (icodcl(ifac,iu).eq.4) then

    ! --- To cancel the mass flux
    isympa(ifac) = 0

    ! Geometric quantities
    srfbnf = surfbn(ifac)

    !===========================================================================
    ! 1. Local framework
    !===========================================================================

    ! Unit normal

    rnx = surfbo(1,ifac)/srfbnf
    rny = surfbo(2,ifac)/srfbnf
    rnz = surfbo(3,ifac)/srfbnf

    ! En ALE, on a eventuellement une vitesse de deplacement de la face
    !   donc seule la composante normale importe (on continue a determiner
    !   TX a partir de la vitesse tangentielle absolue car l'orientation
    !   de TX et T2X est sans importance pour les symetries)
    rcodcn = 0.d0
    if (iale.eq.1) then
      rcodcn = rcodcl(ifac,iu,1)*rnx                           &
             + rcodcl(ifac,iv,1)*rny                           &
             + rcodcl(ifac,iw,1)*rnz
    endif

    upx = velipb(ifac,1)
    upy = velipb(ifac,2)
    upz = velipb(ifac,3)

    if (itytur.eq.3) then

      ! Relative tangential velocity

      usn = upx*rnx+upy*rny+upz*rnz
      tx  = upx -usn*rnx
      ty  = upy -usn*rny
      tz  = upz -usn*rnz
      txn = sqrt( tx**2 +ty**2 +tz**2 )

      ! Unit tangent

      if( txn.ge.epzero) then

        tx = tx/txn
        ty = ty/txn
        tz = tz/txn

      else

        ! If the velocity is zero, vector T is normal and random;
        !   we need it for the reference change for Rij, and we cancel the velocity.

        if(abs(rny).ge.epzero.or.abs(rnz).ge.epzero)then
          rxnn = sqrt(rny**2+rnz**2)
          tx  =  0.d0
          ty  =  rnz/rxnn
          tz  = -rny/rxnn
        elseif(abs(rnx).ge.epzero.or.abs(rnz).ge.epzero)then
          rxnn = sqrt(rnx**2+rnz**2)
          tx  =  rnz/rxnn
          ty  =  0.d0
          tz  = -rnx/rxnn
        else
          write(nfecra,1000)ifac,rnx,rny,rnz
          call csexit (1)
        endif

      endif


      ! --> T2 = RN X T (where X is the cross product)

      t2x = rny*tz - rnz*ty
      t2y = rnz*tx - rnx*tz
      t2z = rnx*ty - rny*tx

      ! --> Orthogonal matrix for change of reference frame ELOGLOij
      !     (from local to global reference frame)

      !                      |TX  -RNX  T2X|
      !             ELOGLO = |TY  -RNY  T2Y|
      !                      |TZ  -RNZ  T2Z|

      !    Its transpose ELOGLOt is its inverse

      eloglo(1,1) =  tx
      eloglo(1,2) = -rnx
      eloglo(1,3) =  t2x
      eloglo(2,1) =  ty
      eloglo(2,2) = -rny
      eloglo(2,3) =  t2y
      eloglo(3,1) =  tz
      eloglo(3,2) = -rnz
      eloglo(3,3) =  t2z

      ! --> Commpute alpha(6,6)

      ! Let f be the center of the boundary faces and
      !   I the center of the matching cell

      ! We noteE Rg (resp. Rl) indexed by f or by I
      !   the Reynolds Stress tensor in the global basis (resp. local)

      ! The alpha matrix applied to the global vector in I'
      !   (Rg11,I'|Rg22,I'|Rg33,I'|Rg12,I'|Rg13,I'|Rg23,I')t
      !    must provide the values to prescribe to the face
      !   (Rg11,f |Rg22,f |Rg33,f |Rg12,f |Rg13,f |Rg23,f )t
      !    except for the Dirichlet boundary conditions (added later)

      ! We define it by computing Rg,f as a function of Rg,I' as follows

      !   RG,f = ELOGLO.RL,f.ELOGLOt (matrix products)

      !                     | RL,I'(1,1)     B*U*.Uk     C*RL,I'(1,3) |
      !      with    RL,f = | B*U*.Uk       RL,I'(2,2)       0        |
      !                     | C*RL,I'(1,3)     0         RL,I'(3,3)   |

      !             with    RL,I = ELOGLOt.RG,I'.ELOGLO
      !                     B = 0
      !              and    C = 0 at the wall (1 with symmetry)

      ! We compute in fact  ELOGLO.projector.ELOGLOt

      clsyme=1.d0
      call clca66 ( clsyme , eloglo , alpha )
      !==========

    endif

    !===========================================================================
    ! 2. Boundary conditions on the velocity
    !    (totaly implicit)
    !    The condition is a zero (except in ALE) Dirichlet on the normal component
    !                     a homogenous Neumann on the other components
    !===========================================================================

    ! Coupled solving of the velocity components

    iel = ifabor(ifac)
    ! --- Physical properties
    visclc = viscl(iel)
    visctc = visct(iel)

    ! --- Geometrical quantity
    distbf = distb(ifac)

    if (itytur.eq.3) then
      hint =   visclc         /distbf
    else
      hint = ( visclc+visctc )/distbf
    endif

    ! Gradient BCs
    coefau(1,ifac) = rcodcn*rnx
    coefau(2,ifac) = rcodcn*rny
    coefau(3,ifac) = rcodcn*rnz

    coefbu(1,1,ifac) = 1.d0-rnx**2
    coefbu(2,2,ifac) = 1.d0-rny**2
    coefbu(3,3,ifac) = 1.d0-rnz**2

    coefbu(1,2,ifac) = -rnx*rny
    coefbu(1,3,ifac) = -rnx*rnz
    coefbu(2,1,ifac) = -rny*rnx
    coefbu(2,3,ifac) = -rny*rnz
    coefbu(3,1,ifac) = -rnz*rnx
    coefbu(3,2,ifac) = -rnz*rny

    ! Flux BCs
    cofafu(1,ifac) = -hint*rcodcn*rnx
    cofafu(2,ifac) = -hint*rcodcn*rny
    cofafu(3,ifac) = -hint*rcodcn*rnz

    cofbfu(1,1,ifac) = hint*rnx**2
    cofbfu(2,2,ifac) = hint*rny**2
    cofbfu(3,3,ifac) = hint*rnz**2

    cofbfu(1,2,ifac) = hint*rnx*rny
    cofbfu(1,3,ifac) = hint*rnx*rnz
    cofbfu(2,1,ifac) = hint*rny*rnx
    cofbfu(2,3,ifac) = hint*rny*rnz
    cofbfu(3,1,ifac) = hint*rnz*rnx
    cofbfu(3,2,ifac) = hint*rnz*rny

    !===========================================================================
    ! 3. Boundary conditions on Rij (partially implicited)
    !===========================================================================

    if (itytur.eq.3) then

      do isou = 1, 6
        fcoefa(isou) = 0.d0
        fcofad(isou) = 0.d0
        fcoefb(isou) = 0.d0
        fcofbd(isou) = 0.d0
      enddo

      do isou = 1,6

        if (isou.eq.1) then
          ivar = ir11
        elseif (isou.eq.2) then
          ivar = ir22
        elseif (isou.eq.3) then
          ivar = ir33
        elseif (isou.eq.4) then
          ivar = ir12
        elseif (isou.eq.5) then
          ivar = ir23
        elseif (isou.eq.6) then
          ivar = ir13
        endif

        ! IMPLICITATION PARTIELLE EVENTUELLE DES CL
        if (iclsyr.eq.1) then
          do ii = 1, 6
            if (ii.ne.isou) then
              fcoefa(isou) = fcoefa(isou) + alpha(isou,ii) * rijipb(ifac,ii)
            endif
          enddo
          fcoefb(isou) = alpha(isou,isou)
        else
          do ii = 1, 6
            fcoefa(isou) = fcoefa(isou) + alpha(isou,ii) * rijipb(ifac,ii)
          enddo
          fcoefb(isou) = 0.d0
        endif

        ! Boundary conditions for the momentum equation
        fcofad(isou) = fcoefa(isou)
        fcofbd(isou) = fcoefb(isou)

        ! Symmetric tensor diffusivity (Daly Harlow -- GGDH)
        if (idften(ivar).eq.6) then

          visci(1,1) = visclc + visten(1,iel)
          visci(2,2) = visclc + visten(2,iel)
          visci(3,3) = visclc + visten(3,iel)
          visci(1,2) =          visten(4,iel)
          visci(2,1) =          visten(4,iel)
          visci(2,3) =          visten(5,iel)
          visci(3,2) =          visten(5,iel)
          visci(1,3) =          visten(6,iel)
          visci(3,1) =          visten(6,iel)

          ! ||Ki.S||^2
          viscis = ( visci(1,1)*surfbo(1,ifac)       &
                   + visci(1,2)*surfbo(2,ifac)       &
                   + visci(1,3)*surfbo(3,ifac))**2   &
                 + ( visci(2,1)*surfbo(1,ifac)       &
                   + visci(2,2)*surfbo(2,ifac)       &
                   + visci(2,3)*surfbo(3,ifac))**2   &
                 + ( visci(3,1)*surfbo(1,ifac)       &
                   + visci(3,2)*surfbo(2,ifac)       &
                   + visci(3,3)*surfbo(3,ifac))**2

          ! IF.Ki.S
          fikis = ( (cdgfbo(1,ifac)-xyzcen(1,iel))*visci(1,1)   &
                  + (cdgfbo(2,ifac)-xyzcen(2,iel))*visci(2,1)   &
                  + (cdgfbo(3,ifac)-xyzcen(3,iel))*visci(3,1)   &
                  )*surfbo(1,ifac)                              &
                + ( (cdgfbo(1,ifac)-xyzcen(1,iel))*visci(1,2)   &
                  + (cdgfbo(2,ifac)-xyzcen(2,iel))*visci(2,2)   &
                  + (cdgfbo(3,ifac)-xyzcen(3,iel))*visci(3,2)   &
                  )*surfbo(2,ifac)                              &
                + ( (cdgfbo(1,ifac)-xyzcen(1,iel))*visci(1,3)   &
                  + (cdgfbo(2,ifac)-xyzcen(2,iel))*visci(2,3)   &
                  + (cdgfbo(3,ifac)-xyzcen(3,iel))*visci(3,3)   &
                  )*surfbo(3,ifac)

          distfi = distb(ifac)

          ! Take I" so that I"F= eps*||FI||*Ki.n when J" is in cell rji
          ! NB: eps =1.d-1 must be consistent with vitens.f90
          fikis = max(fikis, 1.d-1*sqrt(viscis)*distfi)

          hint = viscis/surfbn(ifac)/fikis

        ! Scalar diffusivity
        else
          hint = (visclc+visctc*csrij/cmu)/distbf
        endif

        ! Translate into Diffusive flux BCs
        fcofaf(isou) = -hint*fcoefa(isou)
        fcofbf(isou) = hint*(1.d0-fcoefb(isou))
      enddo
      do isou = 1, 6
        if (irijco.eq.1) then
          coefa_rij(isou, ifac)       = fcoefa(isou)
          coefb_rij(isou,isou, ifac)  = fcoefb(isou)
          coefaf_rij(isou, ifac)      = fcofaf(isou)
          coefbf_rij(isou,isou, ifac) = fcofbf(isou)
          coefad_rij(isou, ifac)      = fcofad(isou)
          coefbd_rij(isou,isou, ifac) = fcofbd(isou)
        else
          if (isou.eq.1) then
            coefa_r11(ifac) = fcoefa(isou)
            coefb_r11(ifac) = fcoefb(isou)
            coefaf_r11(ifac) = fcofaf(isou)
            coefbf_r11(ifac) = fcofbf(isou)
            coefad_r11(ifac) = fcofad(isou)
            coefbd_r11(ifac) = fcofbd(isou)
          else if (isou.eq.2) then
            coefa_r22(ifac) = fcoefa(isou)
            coefb_r22(ifac) = fcoefb(isou)
            coefaf_r22(ifac) = fcofaf(isou)
            coefbf_r22(ifac) = fcofbf(isou)
            coefad_r22(ifac) = fcofad(isou)
            coefbd_r22(ifac) = fcofbd(isou)
          else if (isou.eq.3) then
            coefa_r33(ifac) = fcoefa(isou)
            coefb_r33(ifac) = fcoefb(isou)
            coefaf_r33(ifac) = fcofaf(isou)
            coefbf_r33(ifac) = fcofbf(isou)
            coefad_r33(ifac) = fcofad(isou)
            coefbd_r33(ifac) = fcofbd(isou)
          else if (isou.eq.4) then
            coefa_r12(ifac) = fcoefa(isou)
            coefb_r12(ifac) = fcoefb(isou)
            coefaf_r12(ifac) = fcofaf(isou)
            coefbf_r12(ifac) = fcofbf(isou)
            coefad_r12(ifac) = fcofad(isou)
            coefbd_r12(ifac) = fcofbd(isou)
          else if (isou.eq.5) then
            coefa_r23(ifac) = fcoefa(isou)
            coefb_r23(ifac) = fcoefb(isou)
            coefaf_r23(ifac) = fcofaf(isou)
            coefbf_r23(ifac) = fcofbf(isou)
            coefad_r23(ifac) = fcofad(isou)
            coefbd_r23(ifac) = fcofbd(isou)
          else if (isou.eq.6) then
            coefa_r13(ifac) = fcoefa(isou)
            coefb_r13(ifac) = fcoefb(isou)
            coefaf_r13(ifac) = fcofaf(isou)
            coefbf_r13(ifac) = fcofbf(isou)
            coefad_r13(ifac) = fcofad(isou)
            coefbd_r13(ifac) = fcofbd(isou)
          endif
        endif
      enddo

    endif

  endif
! --- Test sur la presence d'une condition de symetrie vitesse : fin

enddo
! ---  End of loop over boundary faces

!===========================================================================
! 3.bis Boundary conditions on u'T'
!===========================================================================

do iscal = 1, nscal

  if (ityturt(iscal).eq.3) then

    call clsyvt_scalar(iscal, icodcl)
    !=================

  endif

enddo

!===============================================================================
! 4. Symmetry boundary conditions for mesh velocity (ALE module)
!===============================================================================

if (iale.eq.1) then

  call field_get_coefa_v(ivarfl(iuma), claale)
  call field_get_coefb_v(ivarfl(iuma), clbale)
  call field_get_coefaf_v(ivarfl(iuma), cfaale)
  call field_get_coefbf_v(ivarfl(iuma), cfbale)

  call field_get_val_s(iprpfl(ivisma(1)), cpro_visma1)
  call field_get_val_s(iprpfl(ivisma(2)), cpro_visma2)
  call field_get_val_s(iprpfl(ivisma(3)), cpro_visma3)

  do ifac = 1, nfabor
    if (icodcl(ifac,iuma).eq.4) then

      iel = ifabor(ifac)

      ! For a sliding boundary, the normal velocity is enforced to zero
      ! whereas the other components have an Homogenous Neumann
      ! NB: no recontruction in I' here

      ! --- Geometrical quantity
      distbf = distb(ifac)
      srfbnf = surfbn(ifac)

      ! --- Physical properties
      vismsh(1) = cpro_visma1(iel)
      vismsh(2) = cpro_visma2(iel)
      vismsh(3) = cpro_visma3(iel)

      hintv(1) = vismsh(1)/distbf
      hintv(2) = vismsh(2)/distbf
      hintv(3) = vismsh(3)/distbf

      ! Unit normal
      rnx = surfbo(1,ifac)/srfbnf
      rny = surfbo(2,ifac)/srfbnf
      rnz = surfbo(3,ifac)/srfbnf

      ! Coupled solving of the velocity components

      ! Gradient BCs
      claale(1,ifac) = 0.d0
      claale(2,ifac) = 0.d0
      claale(3,ifac) = 0.d0

      clbale(1,1,ifac) = 1.d0-rnx**2
      clbale(2,2,ifac) = 1.d0-rny**2
      clbale(3,3,ifac) = 1.d0-rnz**2

      clbale(1,2,ifac) = -rnx*rny
      clbale(2,1,ifac) = -rny*rnx
      clbale(1,3,ifac) = -rnx*rnz
      clbale(3,1,ifac) = -rnz*rnx
      clbale(2,3,ifac) = -rny*rnz
      clbale(3,2,ifac) = -rnz*rny

      ! Flux BCs
      cfaale(1,ifac) = 0.d0
      cfaale(2,ifac) = 0.d0
      cfaale(3,ifac) = 0.d0

      cfbale(1,1,ifac) = hintv(1)*rnx**2
      cfbale(2,2,ifac) = hintv(2)*rny**2
      cfbale(3,3,ifac) = hintv(3)*rnz**2

      cfbale(1,2,ifac) = hintv(1)*rnx*rny
      cfbale(2,1,ifac) = hintv(2)*rny*rnx
      cfbale(1,3,ifac) = hintv(1)*rnx*rnz
      cfbale(3,1,ifac) = hintv(3)*rnz*rnx
      cfbale(2,3,ifac) = hintv(2)*rny*rnz
      cfbale(3,2,ifac) = hintv(3)*rnz*rny

    endif
  enddo

endif

!===============================================================================
! 5. Formats
!===============================================================================

#if defined(_CS_LANG_FR)

 1000 format(/,' LA NORMALE A LA FACE DE BORD DE SYMETRIE ',I10,/,&
         ' EST NULLE ; COORDONNEES : ',3E12.5)

#else

 1000 format(/,' THE NORMAL TO THE SYMMETRY BOUNDARY FACE ',I10,/,&
         ' IS NULL; COORDINATES: ',3E12.5)

#endif

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

return
end subroutine

!===============================================================================
! Local functions
!===============================================================================

!-------------------------------------------------------------------------------
! Arguments
!______________________________________________________________________________.
!  mode           name          role                                           !
!______________________________________________________________________________!
!> \param[in]     iscal         scalar id
!> \param[in,out] icodcl        face boundary condition code:
!>                               - 1 Dirichlet
!>                               - 3 Neumann
!>                               - 4 sliding and
!>                                 \f$ \vect{u} \cdot \vect{n} = 0 \f$
!>                               - 5 smooth wall and
!>                                 \f$ \vect{u} \cdot \vect{n} = 0 \f$
!>                               - 6 rough wall and
!>                                 \f$ \vect{u} \cdot \vect{n} = 0 \f$
!>                               - 9 free inlet/outlet
!>                                 (input mass flux blocked to 0)
!_______________________________________________________________________________

subroutine clsyvt_scalar &
 ( iscal  , icodcl )

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

use paramx
use numvar
use optcal
use cstphy
use cstnum
use pointe
use entsor
use albase
use parall
use ppppar
use ppthch
use ppincl
use cplsat
use mesh
use field

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

implicit none

! Arguments

integer          iscal

integer          icodcl(nfabor,nvarcl)

! Local variables

integer          ivar, f_id
integer          ifac, iel, isou, jsou
integer          ipccv, ifcvsl

double precision cpp, rkl, visclc
double precision distbf, srfbnf
double precision rnx, rny, rnz
double precision hintt(6)

character(len=80) :: fname

double precision, dimension(:), pointer :: val_s, crom
double precision, dimension(:), pointer :: coefap, coefbp, cofafp, cofbfp
double precision, dimension(:,:), pointer :: coefaut, cofafut, cofarut, visten
double precision, dimension(:,:,:), pointer :: coefbut, cofbfut, cofbrut

double precision, dimension(:), pointer :: viscl, viscls, cpro_cp

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

if (ityturt(iscal).ne.3) return

ivar = isca(iscal)
f_id = ivarfl(ivar)

call field_get_val_s(ivarfl(ivar), val_s)
call field_get_val_v(ivsten, visten)
call field_get_val_s(iprpfl(iviscl), viscl)

call field_get_coefa_s(f_id, coefap)
call field_get_coefb_s(f_id, coefbp)
call field_get_coefaf_s(f_id, cofafp)
call field_get_coefbf_s(f_id, cofbfp)

call field_get_val_s(icrom, crom)

if (icp.gt.0) then
  call field_get_val_s(iprpfl(icp), cpro_cp)
endif

call field_get_key_int (f_id, kivisl, ifcvsl)
if (ifcvsl .ge. 0) then
  call field_get_val_s(ifcvsl, viscls)
endif

ipccv = 0
if (ippmod(icompf) .ge. 0) then
  if (icv.gt.0) then
    ipccv  = ipproc(icv)
  else
    ipccv = 0
  endif
endif

! Turbulent diffusive flux of the scalar T
! (blending factor so that the component v'T' have only
!  mu_T/(mu+mu_T)* Phi_T)

! Name of the scalar ivar
call field_get_name(ivarfl(ivar), fname)

! Index of the corresponding turbulent flux
call field_get_id(trim(fname)//'_turbulent_flux', f_id)

call field_get_coefa_v(f_id,coefaut)
call field_get_coefb_v(f_id,coefbut)
call field_get_coefaf_v(f_id,cofafut)
call field_get_coefbf_v(f_id,cofbfut)
call field_get_coefad_v(f_id,cofarut)
call field_get_coefbd_v(f_id,cofbrut)

! --- Loop on boundary faces
do ifac = 1, nfabor

  ! Test on symmetry boundary condition: start
  if (icodcl(ifac,iu).eq.4) then

    iel = ifabor(ifac)
    ! --- Physical Propreties
    visclc = viscl(iel)
    cpp = 1.d0
    if (iscacp(iscal).eq.1) then
      if (icp.gt.0) then
        cpp = cpro_cp(iel)
      else
        cpp = cp0
      endif
    endif

    ! --- Geometrical quantities
    distbf = distb(ifac)
    srfbnf = surfbn(ifac)

    rnx = surfbo(1,ifac)/srfbnf
    rny = surfbo(2,ifac)/srfbnf
    rnz = surfbo(3,ifac)/srfbnf

    if (ifcvsl .lt. 0) then
      rkl = visls0(iscal)/cpp
    else
      rkl = viscls(iel)/cpp
    endif

    hintt(1) = 0.5d0*(visclc+rkl)/distbf                        &
             + visten(1,iel)*ctheta(iscal)/distbf
    hintt(2) = 0.5d0*(visclc+rkl)/distbf                        &
             + visten(2,iel)*ctheta(iscal)/distbf
    hintt(3) = 0.5d0*(visclc+rkl)/distbf                        &
             + visten(3,iel)*ctheta(iscal)/distbf
    hintt(4) = visten(4,iel)*ctheta(iscal)/distbf
    hintt(5) = visten(5,iel)*ctheta(iscal)/distbf
    hintt(6) = visten(6,iel)*ctheta(iscal)/distbf

    ! Gradient BCs
    coefaut(1,ifac) = 0.d0
    coefaut(2,ifac) = 0.d0
    coefaut(3,ifac) = 0.d0

    coefbut(1,1,ifac) = 1.d0-rnx**2
    coefbut(2,2,ifac) = 1.d0-rny**2
    coefbut(3,3,ifac) = 1.d0-rnz**2

    coefbut(1,2,ifac) = -rnx*rny
    coefbut(1,3,ifac) = -rnx*rnz
    coefbut(2,1,ifac) = -rny*rnx
    coefbut(2,3,ifac) = -rny*rnz
    coefbut(3,1,ifac) = -rnz*rnx
    coefbut(3,2,ifac) = -rnz*rny

    ! Flux BCs
    cofafut(1,ifac) = 0.d0
    cofafut(2,ifac) = 0.d0
    cofafut(3,ifac) = 0.d0

    cofbfut(1,1,ifac) = hintt(1)*rnx**2  + hintt(4)*rnx*rny + hintt(6)*rnx*rnz
    cofbfut(2,2,ifac) = hintt(4)*rnx*rny + hintt(2)*rny**2  + hintt(5)*rny*rnz
    cofbfut(3,3,ifac) = hintt(6)*rnx*rnz + hintt(5)*rny*rnz + hintt(3)*rnz**2

    cofbfut(1,2,ifac) = hintt(1)*rnx*rny + hintt(4)*rny**2  + hintt(6)*rny*rnz
    cofbfut(2,1,ifac) = hintt(1)*rnx*rny + hintt(4)*rny**2  + hintt(6)*rny*rnz
    cofbfut(1,3,ifac) = hintt(1)*rnx*rnz + hintt(4)*rny*rnz + hintt(6)*rnz**2
    cofbfut(3,1,ifac) = hintt(1)*rnx*rnz + hintt(4)*rny*rnz + hintt(6)*rnz**2
    cofbfut(2,3,ifac) = hintt(4)*rnx*rnz + hintt(2)*rny*rnz + hintt(5)*rnz**2
    cofbfut(3,2,ifac) = hintt(4)*rnx*rnz + hintt(2)*rny*rnz + hintt(5)*rnz**2

    ! Boundary conditions for thermal transport equation
    do isou = 1, 3
      cofarut(isou,ifac) = coefaut(isou,ifac)
      do jsou =1, 3
        cofbrut(isou,jsou,ifac) = coefbut(isou,jsou,ifac)
      enddo
    enddo

  endif
  ! Test on velocity symmetry condition: end

enddo

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

return
end subroutine