forked from Nek5000-deprecated/Nek5000-deprecated
-
Notifications
You must be signed in to change notification settings - Fork 2
/
Copy pathbdry.f
2110 lines (2046 loc) · 61.9 KB
/
bdry.f
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
c-----------------------------------------------------------------------
SUBROUTINE SETLOG
C
C Subroutine to initialize logical flags
C
INCLUDE 'SIZE'
INCLUDE 'GEOM'
INCLUDE 'INPUT'
INCLUDE 'TSTEP'
INCLUDE 'TURBO'
INCLUDE 'CTIMER'
COMMON /CPRINT/ IFPRINT
C
common /nekcb/ cb
CHARACTER CB*3
LOGICAL IFALGN,IFNORX,IFNORY,IFNORZ,IFPRINT
C
NFACE = 2*NDIM
NMXV = NFACE*NELV
NMXT = NFACE*NELT
C
IFPRINT = .TRUE.
IFVCOR = .TRUE.
IFGEOM = .FALSE.
IFINTQ = .FALSE.
IFSURT = .FALSE.
IFWCNO = .FALSE.
IFSWALL = .FALSE.
DO 10 IFIELD=1,NFIELD
IFNONL(IFIELD) = .FALSE.
10 CONTINUE
C
CALL LFALSE (IFEPPM,NMXV)
CALL LFALSE (IFQINP,NMXV)
C
IF (IFMODEL) CALL SETSHL
C
IF (IFMVBD) THEN
IFGEOM = .TRUE.
IF ( IFFLOW .AND. .NOT.IFNAV ) IFWCNO = .TRUE.
IF ( IFMELT .AND. .NOT.IFFLOW ) IFWCNO = .TRUE.
ENDIF
C
IF (IFFLOW) THEN
csk call check_cyclic ! fow now; set in .rea file
IFIELD = 1
DO 100 IEL=1,NELV
DO 100 IFC=1,NFACE
CB = CBC(IFC,IEL,IFIELD)
CALL CHKNORD (IFALGN,IFNORX,IFNORY,IFNORZ,IFC,IEL)
IF ( .NOT.IFSTRS ) CALL CHKCBC (CB,IEL,IFC,IFALGN)
IF (CB.EQ.'O ' .OR. CB.EQ.'o ' .OR.
$ CB.EQ.'ON ' .OR. CB.EQ.'on ' .OR.
$ CB.EQ.'S ' .OR. CB.EQ.'s ' .OR.
$ CB.EQ.'SL ' .OR. CB.EQ.'sl ' .OR.
$ CB.EQ.'MM ' .OR. CB.EQ.'mm ' .OR.
$ CB.EQ.'MS ' .OR. CB.EQ.'ms ') THEN
IFVCOR = .FALSE.
IFEPPM(IFC,IEL) = .TRUE.
ENDIF
IF (CB.EQ.'VL ' .OR. CB.EQ.'vl ' .OR.
$ CB.EQ.'WSL' .OR. CB.EQ.'wsl' .OR.
$ CB.EQ.'SL ' .OR. CB.EQ.'sl ' .OR.
$ CB.EQ.'SHL' .OR. CB.EQ.'shl' .OR.
$ CB.EQ.'MM ' .OR. CB.EQ.'mm ' .OR.
$ CB.EQ.'MS ' .OR. CB.EQ.'ms ' .OR.
$ CB.EQ.'O ' .OR. CB.EQ.'o ' .OR.
$ CB.EQ.'ON ' .OR. CB.EQ.'on ') THEN
IFQINP(IFC,IEL) = .TRUE.
ENDIF
IF (CB.EQ.'MS ' .OR. CB.EQ.'ms ' .OR.
$ CB.EQ.'MM ' .OR. CB.EQ.'mm ' .OR.
$ CB.EQ.'MSI' .OR. CB.EQ.'msi' ) THEN
IFSURT = .TRUE.
ENDIF
IF (CB.EQ.'WS ' .OR. CB.EQ.'ws ' .OR.
$ CB.EQ.'WSL' .OR. CB.EQ.'wsl') THEN
IFSWALL = .TRUE.
IFCWUZ = .TRUE.
ENDIF
100 CONTINUE
ENDIF
C
IF (IFHEAT) THEN
C
DO 250 IFIELD=2,NFIELD
DO 250 IEL=1,NELFLD(IFIELD)
DO 250 IFC=1,NFACE
CB=CBC(IFC,IEL,IFIELD)
IF (CB.EQ.'r ' .OR. CB.EQ.'R ') THEN
IFNONL(IFIELD) = .TRUE.
ENDIF
250 CONTINUE
C
ENDIF
if (ifmhd) call set_ifbcor
C
IF (NHIS.GT.0) THEN
IQ = 0
DO 300 IH=1,NHIS
IF ( HCODE(10,IH) .EQ. 'I' ) THEN
IFINTQ = .TRUE.
IOBJ = LOCHIS(1,IH)
IQ = IQ + 1
IF (IOBJ.GT.NOBJ .OR. IOBJ.LT.0) THEN
WRITE (6,*)
$ 'ERROR : Undefined Object for integral',IQ
call exitt
ENDIF
ENDIF
300 CONTINUE
ENDIF
C
C Establish global consistency of LOGICALS amongst all processors.
C
CALL GLLOG(IFVCOR , .FALSE.)
CALL GLLOG(IFSURT , .TRUE. )
CALL GLLOG(IFSWALL, .TRUE. )
CALL GLLOG(IFCWUZ , .TRUE. )
CALL GLLOG(IFWCNO , .TRUE. )
DO 400 IFIELD=2,NFIELD
CALL GLLOG(IFNONL(IFIELD),.TRUE.)
400 CONTINUE
C
IF (NIO.EQ.0) THEN
WRITE (6,*) 'IFTRAN =',IFTRAN
WRITE (6,*) 'IFFLOW =',IFFLOW
WRITE (6,*) 'IFHEAT =',IFHEAT
WRITE (6,*) 'IFSPLIT =',IFSPLIT
WRITE (6,*) 'IFLOMACH =',IFLOMACH
WRITE (6,*) 'IFUSERVP =',IFUSERVP
WRITE (6,*) 'IFUSERMV =',IFUSERMV
WRITE (6,*) 'IFSTRS =',IFSTRS
WRITE (6,*) 'IFCHAR =',IFCHAR
WRITE (6,*) 'IFCYCLIC =',IFCYCLIC
WRITE (6,*) 'IFAXIS =',IFAXIS
WRITE (6,*) 'IFMVBD =',IFMVBD
WRITE (6,*) 'IFMELT =',IFMELT
WRITE (6,*) 'IFMODEL =',IFMODEL
WRITE (6,*) 'IFKEPS =',IFKEPS
WRITE (6,*) 'IFMOAB =',IFMOAB
WRITE (6,*) 'IFNEKNEK =',IFNEKNEK
WRITE (6,*) 'IFSYNC =',IFSYNC
WRITE (6,*) ' '
WRITE (6,*) 'IFVCOR =',IFVCOR
WRITE (6,*) 'IFINTQ =',IFINTQ
WRITE (6,*) 'IFCWUZ =',IFCWUZ
WRITE (6,*) 'IFSWALL =',IFSWALL
WRITE (6,*) 'IFGEOM =',IFGEOM
WRITE (6,*) 'IFSURT =',IFSURT
WRITE (6,*) 'IFWCNO =',IFWCNO
WRITE (6,*) 'IFCMT =',IFCMT
WRITE (6,*) 'IFVISC =',IFVISC
WRITE (6,*) 'IFFLTR =',IFFLTR
DO 500 IFIELD=1,NFIELD
WRITE (6,*) ' '
WRITE (6,*) 'IFTMSH for field',IFIELD,' = ',IFTMSH(IFIELD)
WRITE (6,*) 'IFADVC for field',IFIELD,' = ',IFADVC(IFIELD)
WRITE (6,*) 'IFNONL for field',IFIELD,' = ',IFNONL(IFIELD)
500 CONTINUE
WRITE (6,*) ' '
if (param(99).gt.0) write(6,*) 'Dealiasing enabled, lxd=', lxd
ENDIF
C
RETURN
END
C
c-----------------------------------------------------------------------
SUBROUTINE SETRZER
C-------------------------------------------------------------------
C
C Check for axisymmetric case.
C Are some of the elements close to the axis?
C
C-------------------------------------------------------------------
INCLUDE 'SIZE'
INCLUDE 'GEOM'
INCLUDE 'INPUT'
C
C Single or double precision???
C
DELTA = 1.E-9
X = 1.+DELTA
Y = 1.
DIFF = ABS(X-Y)
IF (DIFF.EQ.0.) EPS = 1.E-7
IF (DIFF.GT.0.) EPS = 1.E-14
eps1 = 1.e-6 ! for prenek mesh in real*4
C
DO 100 IEL=1,NELT
IFRZER(IEL) = .FALSE.
IF (IFAXIS) THEN
NVERT = 0
DO 10 IC=1,4
IF(ABS(YC(IC,IEL)).LT.EPS1) THEN
NVERT = NVERT+1
YC(IC,IEL) = 0.0 ! exactly on the axis
ENDIF
10 CONTINUE
ENDIF
IEDGE = 1
IF ((NVERT.EQ.2).AND.(CCURVE(IEDGE,IEL).EQ.' '))
$ IFRZER(IEL) = .TRUE.
100 CONTINUE
RETURN
END
C
c-----------------------------------------------------------------------
SUBROUTINE CHKNORD (IFALGN,IFNORX,IFNORY,IFNORZ,IFC,IEL)
C
C Check direction of normal of an element face for
C alignment with the X, Y, or Z axis.
C
INCLUDE 'SIZE'
INCLUDE 'GEOM'
C
LOGICAL IFALGN,IFNORX,IFNORY,IFNORZ
C
SUMX = 0.0
SUMY = 0.0
SUMZ = 0.0
TOLNOR = 1.0e-3
IFALGN = .FALSE.
IFNORX = .FALSE.
IFNORY = .FALSE.
IFNORZ = .FALSE.
C
IF (NDIM.EQ.2) THEN
C
NCPF = NX1
DO 100 IX=1,NX1
SUMX = SUMX + ABS( ABS(UNX(IX,1,IFC,IEL)) - 1.0 )
SUMY = SUMY + ABS( ABS(UNY(IX,1,IFC,IEL)) - 1.0 )
100 CONTINUE
SUMX = SUMX / NCPF
SUMY = SUMY / NCPF
IF ( SUMX.LT.TOLNOR ) THEN
IFNORX = .TRUE.
IFALGN = .TRUE.
ENDIF
IF ( SUMY.LT.TOLNOR ) THEN
IFNORY = .TRUE.
IFALGN = .TRUE.
ENDIF
C
ELSE
C
NCPF = NX1*NX1
DO 200 IX=1,NX1
DO 200 IY=1,NY1
SUMX = SUMX + ABS( ABS(UNX(IX,IY,IFC,IEL)) - 1.0 )
SUMY = SUMY + ABS( ABS(UNY(IX,IY,IFC,IEL)) - 1.0 )
SUMZ = SUMZ + ABS( ABS(UNZ(IX,IY,IFC,IEL)) - 1.0 )
200 CONTINUE
SUMX = SUMX / NCPF
SUMY = SUMY / NCPF
SUMZ = SUMZ / NCPF
IF ( SUMX.LT.TOLNOR ) THEN
IFNORX = .TRUE.
IFALGN = .TRUE.
ENDIF
IF ( SUMY.LT.TOLNOR ) THEN
IFNORY = .TRUE.
IFALGN = .TRUE.
ENDIF
IF ( SUMZ.LT.TOLNOR ) THEN
IFNORZ = .TRUE.
IFALGN = .TRUE.
ENDIF
C
ENDIF
C
RETURN
END
c-----------------------------------------------------------------------
SUBROUTINE CHKAXCB
C
INCLUDE 'SIZE'
INCLUDE 'INPUT'
CHARACTER CB*3
C
IFLD = 1
NFACE = 2*NDIM
C
DO 100 IEL=1,NELV
DO 100 IFC=1,NFACE
CB = CBC(IFC,IEL,IFLD)
IF (CB.EQ.'A ' .AND. IFC.NE.1) GOTO 9000
100 CONTINUE
C
RETURN
C
9000 WRITE (6,*) ' Element face on the axis of symmetry must be FACE 1'
WRITE (6,*) ' Element',IEL,' face',IFC,' is on the axis.'
call exitt
C
END
c-----------------------------------------------------------------------
SUBROUTINE CHKCBC (CB,IEL,IFC,IFALGN)
include 'SIZE'
include 'PARALLEL'
C
C Check for illegal boundary conditions
C
CHARACTER CB*3
LOGICAL IFALGN
ieg = lglel(iel)
C
C Laplacian formulation only
C
IF (CB.EQ.'SH ' .OR. CB.EQ.'sh ' .OR.
$ CB.EQ.'SHL' .OR. CB.EQ.'shl' .OR.
$ CB.EQ.'S ' .OR. CB.EQ.'s ' .OR.
$ CB.EQ.'SL ' .OR. CB.EQ.'sl ' .OR.
$ CB.EQ.'MM ' .OR. CB.EQ.'mm ' .OR.
$ CB.EQ.'MS ' .OR. CB.EQ.'ms ' .OR.
$ CB.EQ.'MSI' .OR. CB.EQ.'msi' ) GOTO 9001
IF ( .NOT.IFALGN .AND.
$ (CB.EQ.'ON ' .OR. CB.EQ.'on ' .OR. CB.EQ.'SYM') ) GOTO 9010
RETURN
C
9001 WRITE (6,*) ' Illegal traction boundary conditions detected for'
GOTO 9999
9010 WRITE (6,*) ' Mixed B.C. on a side nonaligned with either the X,Y,
$ or Z axis detected for'
9999 WRITE (6,*) ' Element',ieg,' side',IFC,'.'
WRITE (6,*) ' Selected option only allowed for STRESS FORMULATION'
WRITE (6,*) ' Execution terminates'
call exitt
END
c-----------------------------------------------------------------------
SUBROUTINE BCMASK
C
C Zero out masks corresponding to Dirichlet boundary points.
C
INCLUDE 'SIZE'
INCLUDE 'TSTEP'
INCLUDE 'INPUT'
INCLUDE 'MVGEOM'
INCLUDE 'SOLN'
INCLUDE 'TOPOL'
common /nekcb/ cb
character*3 cb
character*1 cb1(3)
equivalence (cb1,cb)
logical ifalgn,ifnorx,ifnory,ifnorz
integer e,f
NFACES=2*NDIM
NXYZ =NX1*NY1*NZ1
C
C Masks for moving mesh
C
IF (IFMVBD) THEN
IFIELD = 0
CALL STSMASK (W1MASK,W2MASK,W3MASK)
do e=1,nelv
do f=1,nfaces
if (cbc(f,e,1).eq.'msi'.or.cbc(f,e,1).eq.'msi') then
call facev(w1mask,e,f,0.0,nx1,ny1,nz1)
call facev(w2mask,e,f,0.0,nx1,ny1,nz1)
call facev(w3mask,e,f,0.0,nx1,ny1,nz1)
endif
enddo
enddo
ENDIF
C
C Masks for flow variables
C
IF (IFFLOW) THEN
IFIELD = 1
NEL = NELFLD(IFIELD)
NTOT = NXYZ*NEL
C
C Pressure mask
C
call rone(pmask,ntot)
do 50 iel=1,nelt
do 50 iface=1,nfaces
cb=cbc(iface,iel,ifield)
if (cb.eq.'O ' .or. cb.eq.'ON ' .or.
$ cb.eq.'o ' .or. cb.eq.'on ')
$ call facev(pmask,iel,iface,0.0,nx1,ny1,nz1)
50 continue
if (nelt.gt.nelv) then
nn=nx1*ny1*nz1*(nelt-nelv)
call rzero(pmask(1,1,1,nelv+1),nn)
endif
C
C Zero out mask at Neumann-Dirichlet interfaces
C
CALL DSOP(PMASK,'MUL',NX1,NY1,NZ1)
C
C Velocity masks
C
c write(6,*) 'MASK ifstrs',ifstrs,ifield
c call exitt
IF (IFSTRS) THEN
CALL STSMASK (V1MASK,V2MASK,V3MASK)
ELSE
C
CALL RONE(V1MASK,NTOT)
CALL RONE(V2MASK,NTOT)
CALL RONE(V3MASK,NTOT)
CALL RONE( OMASK,NTOT)
C
DO 100 IEL=1,NELV
DO 100 IFACE=1,NFACES
CB =CBC(IFACE,IEL,IFIELD)
CALL CHKNORD (IFALGN,IFNORX,IFNORY,IFNORZ,IFACE,IEL)
C
C All-Dirichlet boundary conditions
C
IF (CB.EQ.'v ' .OR. CB.EQ.'V ' .OR. CB.EQ.'vl ' .OR.
$ CB.EQ.'VL ' .OR. CB.EQ.'W ') THEN
CALL FACEV (V1MASK,IEL,IFACE,0.0,NX1,NY1,NZ1)
CALL FACEV (V2MASK,IEL,IFACE,0.0,NX1,NY1,NZ1)
CALL FACEV (V3MASK,IEL,IFACE,0.0,NX1,NY1,NZ1)
GOTO 100
ENDIF
C
C Mixed-Dirichlet-Neumann boundary conditions
C
IF (CB.EQ.'SYM') THEN
IF ( .NOT.IFALGN .OR. IFNORX )
$ CALL FACEV (V1MASK,IEL,IFACE,0.0,NX1,NY1,NZ1)
IF ( IFNORY )
$ CALL FACEV (V2MASK,IEL,IFACE,0.0,NX1,NY1,NZ1)
IF ( IFNORZ )
$ CALL FACEV (V3MASK,IEL,IFACE,0.0,NX1,NY1,NZ1)
GOTO 100
ENDIF
IF (CB.EQ.'ON ') THEN
IF ( IFNORY .OR. IFNORZ )
$ CALL FACEV (V1MASK,IEL,IFACE,0.0,NX1,NY1,NZ1)
IF ( .NOT.IFALGN .OR. IFNORX .OR. IFNORZ )
$ CALL FACEV (V2MASK,IEL,IFACE,0.0,NX1,NY1,NZ1)
IF ( .NOT.IFALGN .OR. IFNORX .OR. IFNORY )
$ CALL FACEV (V3MASK,IEL,IFACE,0.0,NX1,NY1,NZ1)
GOTO 100
ENDIF
IF (CB.EQ.'A ') THEN
CALL FACEV (V2MASK,IEL,IFACE,0.0,NX1,NY1,NZ1)
CALL FACEV ( OMASK,IEL,IFACE,0.0,NX1,NY1,NZ1)
ENDIF
100 CONTINUE
CALL DSOP ( OMASK,'MUL',NX1,NY1,NZ1)
call opdsop(v1mask,v2mask,v3mask,'MUL') ! no rotation for mul
ENDIF
C
ENDIF
C
C Masks for passive scalars +
C k and e if k-e turbulence modem:
C k = nfield-1
C e = nfield
C
IF (IFHEAT) THEN
C
DO 1200 IFIELD=2,NFIELD
IPSCAL = IFIELD-1
NEL = NELFLD(IFIELD)
NTOT = NXYZ*NEL
CALL RONE (TMASK(1,1,1,1,IPSCAL),NTOT)
DO 1100 IEL=1,NEL
DO 1100 IFACE=1,NFACES
CB =CBC(IFACE,IEL,IFIELD)
C
C Assign mask values.
C
IF (CB.EQ.'T ' .OR. CB.EQ.'t ' .OR.
$ (CB.EQ.'A ' .AND. IFAZIV) .OR.
$ CB.EQ.'MCI' .OR. CB.EQ.'MLI' .OR.
$ CB.EQ.'KD ' .OR. CB.EQ.'kd ' .OR.
$ CB.EQ.'ED ' .OR. CB.EQ.'ed ' .OR.
$ CB.EQ.'KW ' .OR. CB.EQ.'KWS' .OR. CB.EQ.'EWS')
$ CALL FACEV (TMASK(1,1,1,1,IPSCAL),
$ IEL,IFACE,0.0,NX1,NY1,NZ1)
1100 CONTINUE
CALL DSOP (TMASK(1,1,1,1,IPSCAL),'MUL',NX1,NY1,NZ1)
1200 CONTINUE
C
ENDIF
C
C Masks for B-field
C
if (ifmhd) then
ifield = ifldmhd
nel = nelfld(ifield)
ntot = nxyz*nel
C
C B-field pressure mask
C
call rone(bpmask,ntot)
do iel=1,nelv
do iface=1,nfaces
cb=cbc(iface,iel,ifield)
if (cb.eq.'O ' .or. cb.eq.'ON ')
$ call facev(bpmask,iel,iface,0.0,nx1,ny1,nz1)
enddo
enddo
C
C Zero out mask at Neumann-Dirichlet interfaces
C
call dsop(bpmask,'MUL',nx1,ny1,nz1)
C
C B-field masks
C
if (ifstrs) then
call stsmask (b1mask,b2mask,b3mask)
else
C
call rone(b1mask,ntot)
call rone(b2mask,ntot)
call rone(b3mask,ntot)
C
do iel=1,nelv
do iface=1,nfaces
cb =cbc(iface,iel,ifield)
call chknord (ifalgn,ifnorx,ifnory,ifnorz,iface,iel)
c
if (cb.eq.'v ' .or. cb.eq.'V ' .or. cb.eq.'vl ' .or.
$ cb.eq.'VL ' .or. cb.eq.'W ') then
c
c All-Dirichlet boundary conditions
c
call facev (b1mask,iel,iface,0.0,nx1,ny1,nz1)
call facev (b2mask,iel,iface,0.0,nx1,ny1,nz1)
call facev (b3mask,iel,iface,0.0,nx1,ny1,nz1)
c
elseif (cb.eq.'SYM') then
c
c Mixed-Dirichlet-Neumann boundary conditions
c
if ( .not.ifalgn .or. ifnorx )
$ call facev (b1mask,iel,iface,0.0,nx1,ny1,nz1)
if ( ifnory )
$ call facev (b2mask,iel,iface,0.0,nx1,ny1,nz1)
if ( ifnorz )
$ call facev (b3mask,iel,iface,0.0,nx1,ny1,nz1)
c
elseif (cb.eq.'ON ') then
c
c Mixed-Dirichlet-Neumann boundary conditions
c
if ( ifnory .or. ifnorz )
$ call facev (b1mask,iel,iface,0.0,nx1,ny1,nz1)
if ( .not.ifalgn .or. ifnorx .or. ifnorz )
$ call facev (b2mask,iel,iface,0.0,nx1,ny1,nz1)
if ( .not.ifalgn .or. ifnorx .or. ifnory )
$ call facev (b3mask,iel,iface,0.0,nx1,ny1,nz1)
c
elseif (cb.eq.'A ') then
c
c axisymmetric centerline
c
call facev (b2mask,iel,iface,0.0,nx1,ny1,nz1)
c
else
c
if ( cb1(1).eq.'d' )
$ call facev (b1mask,iel,iface,0.0,nx1,ny1,nz1)
if ( cb1(2).eq.'d' )
$ call facev (b2mask,iel,iface,0.0,nx1,ny1,nz1)
if ( cb1(3).eq.'d' .and. if3d )
$ call facev (b3mask,iel,iface,0.0,nx1,ny1,nz1)
c
endif
enddo
enddo
c
call dsop(b1mask,'MUL',nx1,ny1,nz1)
call dsop(b2mask,'MUL',nx1,ny1,nz1)
if (ndim.eq.3) call dsop(b3mask,'MUL',nx1,ny1,nz1)
endif
endif
C
RETURN
END
c-----------------------------------------------------------------------
SUBROUTINE BCDIRVC(V1,V2,V3,mask1,mask2,mask3)
C
C Apply Dirichlet boundary conditions to surface of vector (V1,V2,V3).
C Use IFIELD as a guide to which boundary conditions are to be applied.
C
INCLUDE 'SIZE'
INCLUDE 'TSTEP'
INCLUDE 'INPUT'
INCLUDE 'GEOM'
INCLUDE 'SOLN'
INCLUDE 'TOPOL'
INCLUDE 'CTIMER'
COMMON /SCRUZ/ TMP1(LX1,LY1,LZ1,LELV)
$ , TMP2(LX1,LY1,LZ1,LELV)
$ , TMP3(LX1,LY1,LZ1,LELV)
COMMON /SCRMG/ TMQ1(LX1,LY1,LZ1,LELV)
$ , TMQ2(LX1,LY1,LZ1,LELV)
$ , TMQ3(LX1,LY1,LZ1,LELV)
C
REAL V1(NX1,NY1,NZ1,LELV),V2(NX1,NY1,NZ1,LELV)
$ ,V3(NX1,NY1,NZ1,LELV)
real mask1(nx1,ny1,nz1,lelv),mask2(nx1,ny1,nz1,lelv)
$ ,mask3(nx1,ny1,nz1,lelv)
c
common /nekcb/ cb
character cb*3
character*1 cb1(3)
equivalence (cb1,cb)
c
logical ifonbc
c
ifonbc = .false.
c
if (icalld.eq.0) then
tusbc=0.0
nusbc=0
icalld=icalld+1
endif
nusbc=nusbc+1
etime1=dnekclock()
C
C
NFACES=2*NDIM
NXYZ =NX1*NY1*NZ1
NEL =NELFLD(IFIELD)
NTOT =NXYZ*NEL
C
CALL RZERO(TMP1,NTOT)
CALL RZERO(TMP2,NTOT)
IF (IF3D) CALL RZERO(TMP3,NTOT)
C
C Velocity boundary conditions
C
c write(6,*) 'BCDIRV: ifield',ifield
DO 2100 ISWEEP=1,2
DO 2000 IE=1,NEL
DO 2000 IFACE=1,NFACES
CB = CBC(IFACE,IE,IFIELD)
BC1 = BC(1,IFACE,IE,IFIELD)
BC2 = BC(2,IFACE,IE,IFIELD)
BC3 = BC(3,IFACE,IE,IFIELD)
IF (CB.EQ.'V ' .OR. CB.EQ.'VL ' .OR.
$ CB.EQ.'WS ' .OR. CB.EQ.'WSL') THEN
CALL FACEV (TMP1,IE,IFACE,BC1,NX1,NY1,NZ1)
CALL FACEV (TMP2,IE,IFACE,BC2,NX1,NY1,NZ1)
IF (IF3D) CALL FACEV (TMP3,IE,IFACE,BC3,NX1,NY1,NZ1)
IF ( IFQINP(IFACE,IE) )
$ CALL GLOBROT (TMP1(1,1,1,IE),TMP2(1,1,1,IE),
$ TMP3(1,1,1,IE),IE,IFACE)
ENDIF
IF (CB.EQ.'v ' .OR. CB.EQ.'vl ' .OR.
$ CB.EQ.'ws ' .OR. CB.EQ.'wsl' .OR.
$ CB.EQ.'mv ' .OR. CB.EQ.'mvn' .OR.
$ cb1(1).eq.'d'.or.cb1(2).eq.'d'.or.cb1(3).eq.'d') then
call faceiv (cb,tmp1(1,1,1,ie),tmp2(1,1,1,ie),
$ tmp3(1,1,1,ie),ie,iface,nx1,ny1,nz1)
IF ( IFQINP(IFACE,IE) )
$ CALL GLOBROT (TMP1(1,1,1,IE),TMP2(1,1,1,IE),
$ TMP3(1,1,1,IE),IE,IFACE)
ENDIF
IF (CB.EQ.'ON ' .OR. CB.EQ.'on ') then ! 5/21/01 pff
ifonbc =.true.
CALL FACEIV ('v ',TMP1(1,1,1,IE),TMP2(1,1,1,IE),
$ TMP3(1,1,1,IE),IE,IFACE,NX1,NY1,NZ1)
ENDIF
2000 CONTINUE
DO 2010 IE=1,NEL
DO 2010 IFACE=1,NFACES
IF (CBC(IFACE,IE,IFIELD).EQ.'W ') THEN
CALL FACEV (TMP1,IE,IFACE,0.0,NX1,NY1,NZ1)
CALL FACEV (TMP2,IE,IFACE,0.0,NX1,NY1,NZ1)
IF (IF3D) CALL FACEV (TMP3,IE,IFACE,0.0,NX1,NY1,NZ1)
ENDIF
2010 CONTINUE
C
C Take care of Neumann-Dirichlet shared edges...
C
if (isweep.eq.1) then
call opdsop(tmp1,tmp2,tmp3,'MXA')
else
call opdsop(tmp1,tmp2,tmp3,'MNA')
endif
2100 CONTINUE
C
C Copy temporary array to velocity arrays.
C
IF ( .NOT.IFSTRS ) THEN
CALL COL2(V1,mask1,NTOT)
CALL COL2(V2,mask2,NTOT)
IF (IF3D) CALL COL2(V3,mask3,NTOT)
if (ifonbc) then
call antimsk1(tmp1,mask1,ntot)
call antimsk1(tmp2,mask2,ntot)
if (if3d) call antimsk1(tmp3,mask3,ntot)
endif
ELSE
IF (IFMODEL) THEN
CALL COPY (TMQ1,TMP1,NTOT)
CALL COPY (TMQ2,TMP2,NTOT)
IF (NDIM.EQ.3) CALL COPY (TMQ3,TMP3,NTOT)
CALL AMASK (TMP1,TMP2,TMP3,TMQ1,TMQ2,TMQ3,NELV)
ENDIF
CALL RMASK (V1,V2,V3,NELV)
ENDIF
CALL ADD2(V1,TMP1,NTOT)
CALL ADD2(V2,TMP2,NTOT)
IF (IF3D) CALL ADD2(V3,TMP3,NTOT)
tusbc=tusbc+(dnekclock()-etime1)
RETURN
END
c-----------------------------------------------------------------------
SUBROUTINE BCDIRSC(S)
C
C Apply Dirichlet boundary conditions to surface of scalar, S.
C Use IFIELD as a guide to which boundary conditions are to be applied.
C
INCLUDE 'SIZE'
INCLUDE 'TSTEP'
INCLUDE 'INPUT'
INCLUDE 'SOLN'
INCLUDE 'TOPOL'
INCLUDE 'CTIMER'
C
DIMENSION S(LX1,LY1,LZ1,LELT)
COMMON /SCRSF/ TMP(LX1,LY1,LZ1,LELT)
$ , TMA(LX1,LY1,LZ1,LELT)
$ , SMU(LX1,LY1,LZ1,LELT)
common /nekcb/ cb
CHARACTER CB*3
if (icalld.eq.0) then
tusbc=0.0
nusbc=0
icalld=icalld+1
endif
nusbc=nusbc+1
etime1=dnekclock()
C
IFLD = 1
NFACES = 2*NDIM
NXYZ = NX1*NY1*NZ1
NEL = NELFLD(IFIELD)
NTOT = NXYZ*NEL
NFLDT = NFIELD - 1
C
CALL RZERO(TMP,NTOT)
C
C Temperature boundary condition
C
DO 2100 ISWEEP=1,2
C
IF (IFMODEL .AND. IFKEPS .AND. IFIELD.GE.NFLDT)
$ CALL TURBWBC (TMP,TMA,SMU)
C
DO 2010 IE=1,NEL
DO 2010 IFACE=1,NFACES
CB=CBC(IFACE,IE,IFIELD)
BC1=BC(1,IFACE,IE,IFIELD)
BC2=BC(2,IFACE,IE,IFIELD)
BC3=BC(3,IFACE,IE,IFIELD)
BC4=BC(4,IFACE,IE,IFIELD)
BCK=BC(4,IFACE,IE,IFLD)
BCE=BC(5,IFACE,IE,IFLD)
IF (CB.EQ.'T ') CALL FACEV (TMP,IE,IFACE,BC1,NX1,NY1,NZ1)
IF (CB.EQ.'MCI') CALL FACEV (TMP,IE,IFACE,BC4,NX1,NY1,NZ1)
IF (CB.EQ.'MLI') CALL FACEV (TMP,IE,IFACE,BC4,NX1,NY1,NZ1)
IF (CB.EQ.'KD ') CALL FACEV (TMP,IE,IFACE,BCK,NX1,NY1,NZ1)
IF (CB.EQ.'ED ') CALL FACEV (TMP,IE,IFACE,BCE,NX1,NY1,NZ1)
IF (CB.EQ.'t ' .OR. CB.EQ.'kd ' .or.
$ CB.EQ.'ed ' .or. cb.eq.'o ')
$ CALL FACEIS (CB,TMP(1,1,1,IE),IE,IFACE,NX1,NY1,NZ1)
2010 CONTINUE
C
C Take care of Neumann-Dirichlet shared edges...
C
IF (ISWEEP.EQ.1) CALL DSOP(TMP,'MXA',NX1,NY1,NZ1)
IF (ISWEEP.EQ.2) CALL DSOP(TMP,'MNA',NX1,NY1,NZ1)
2100 CONTINUE
C
C Copy temporary array to temperature array.
C
CALL COL2(S,TMASK(1,1,1,1,IFIELD-1),NTOT)
CALL ADD2(S,TMP,NTOT)
tusbc=tusbc+(dnekclock()-etime1)
RETURN
END
C
c-----------------------------------------------------------------------
SUBROUTINE BCNEUSC(S,ITYPE)
C
C Apply Neumann boundary conditions to surface of scalar, S.
C Use IFIELD as a guide to which boundary conditions are to be applied.
C
C If ITYPE = 1, then S is returned as the rhs contribution to the
C volumetric flux.
C
C If ITYPE =-1, then S is returned as the lhs contribution to the
C diagonal of A.
C
C
INCLUDE 'SIZE'
INCLUDE 'TOTAL'
INCLUDE 'CTIMER'
INCLUDE 'NEKUSE'
C
DIMENSION S(LX1,LY1,LZ1,LELT)
common /nekcb/ cb
CHARACTER CB*3
C
if (icalld.eq.0) then
tusbc=0.0
nusbc=0
icalld=icalld+1
endif
nusbc=nusbc+1
etime1=dnekclock()
C
NFACES=2*NDIM
NXYZ =NX1*NY1*NZ1
NEL =NELFLD(IFIELD)
NTOT =NXYZ*NEL
CALL RZERO(S,NTOT)
C
IF (ITYPE.EQ.-1) THEN
C
C Compute diagonal contributions to accomodate Robin boundary conditions
C
DO 1000 IE=1,NEL
DO 1000 IFACE=1,NFACES
ieg=lglel(ie)
CB =CBC(IFACE,IE,IFIELD)
IF (CB.EQ.'C ' .OR. CB.EQ.'c ' .OR.
$ CB.EQ.'R ' .OR. CB.EQ.'r ') THEN
C
IF (CB.EQ.'C ') HC = BC(2,IFACE,IE,IFIELD)
IF (CB.EQ.'R ') THEN
TINF = BC(1,IFACE,IE,IFIELD)
HRAD = BC(2,IFACE,IE,IFIELD)
ENDIF
IA=0
C
C IA is areal counter, assumes advancing fastest index first. (IX...IY...IZ)
C
CALL FACIND (KX1,KX2,KY1,KY2,KZ1,KZ2,NX1,NY1,NZ1,IFACE)
DO 100 IZ=KZ1,KZ2
DO 100 IY=KY1,KY2
DO 100 IX=KX1,KX2
IA = IA + 1
TS = T(IX,IY,IZ,IE,IFIELD-1)
IF (CB.EQ.'c ' .OR. CB.EQ.'r ') THEN
CALL NEKASGN (IX,IY,IZ,IE)
CALL USERBC (IX,IY,IZ,IFACE,IEG)
ENDIF
IF (CB.EQ.'r ' .OR. CB.EQ.'R ')
$ HC = HRAD * (TINF**2 + TS**2) * (TINF + TS)
S(IX,IY,IZ,IE) = S(IX,IY,IZ,IE) +
$ HC*AREA(IA,1,IFACE,IE)/BM1(IX,IY,IZ,IE)
100 CONTINUE
ENDIF
1000 CONTINUE
ENDIF
IF (ITYPE.EQ.1) THEN
C
C Add passive scalar fluxes to rhs
C
DO 2000 IE=1,NEL
DO 2000 IFACE=1,NFACES
ieg=lglel(ie)
CB =CBC(IFACE,IE,IFIELD)
IF (CB.EQ.'F ' .OR. CB.EQ.'f ' .OR.
$ CB.EQ.'C ' .OR. CB.EQ.'c ' .OR.
$ CB.EQ.'R ' .OR. CB.EQ.'r ' ) THEN
C
IF (CB.EQ.'F ') FLUX=BC(1,IFACE,IE,IFIELD)
IF (CB.EQ.'C ') FLUX=BC(1,IFACE,IE,IFIELD)
$ *BC(2,IFACE,IE,IFIELD)
IF (CB.EQ.'R ') THEN
TINF=BC(1,IFACE,IE,IFIELD)
HRAD=BC(2,IFACE,IE,IFIELD)
ENDIF
C
C Add local weighted flux values to rhs, S.
C
C IA is areal counter, assumes advancing fastest index first. (IX...IY...IZ)
IA=0
CALL FACIND (KX1,KX2,KY1,KY2,KZ1,KZ2,NX1,NY1,NZ1,IFACE)
DO 200 IZ=KZ1,KZ2
DO 200 IY=KY1,KY2
DO 200 IX=KX1,KX2
IA = IA + 1
TS = T(IX,IY,IZ,IE,IFIELD-1)
IF (CB.EQ.'f ') THEN
CALL NEKASGN (IX,IY,IZ,IE)
CALL USERBC (IX,IY,IZ,IFACE,IEG)
ENDIF
IF (CB.EQ.'c ') THEN
CALL NEKASGN (IX,IY,IZ,IE)
CALL USERBC (IX,IY,IZ,IFACE,IEG)
FLUX = TINF*HC
ENDIF
IF (CB.EQ.'r ') THEN
CALL NEKASGN (IX,IY,IZ,IE)
CALL USERBC (IX,IY,IZ,IFACE,IEG)
ENDIF
IF (CB.EQ.'R ' .OR. CB.EQ.'r ')
$ FLUX = HRAD*(TINF**2 + TS**2)*(TINF + TS) * TINF
C
C Add computed fluxes to boundary surfaces:
C
S(IX,IY,IZ,IE) = S(IX,IY,IZ,IE)
$ + FLUX*AREA(IA,1,IFACE,IE)
200 CONTINUE
ENDIF
2000 CONTINUE
ENDIF
C
tusbc=tusbc+(dnekclock()-etime1)
C
RETURN
END
c-----------------------------------------------------------------------
SUBROUTINE FACEIS (CB,S,IEL,IFACE,NX,NY,NZ)
C
C Assign inflow boundary conditions to face(IE,IFACE)
C for scalar S.
C
INCLUDE 'SIZE'
INCLUDE 'PARALLEL'
INCLUDE 'NEKUSE'
INCLUDE 'TSTEP' ! ifield 11/19/2010
INCLUDE 'SOLN' ! tmask() 11/19/2010
C
DIMENSION S(LX1,LY1,LZ1)
CHARACTER CB*3
c
common /nekcb/ cb3
character*3 cb3
cb3 = cb
ifld1 = ifield-1
C Passive scalar term
ieg = lglel(iel)
CALL FACIND (KX1,KX2,KY1,KY2,KZ1,KZ2,NX,NY,NZ,IFACE)
if (cb.eq.'t ') then
DO 100 IZ=KZ1,KZ2 ! 11/19/2010: The tmask() screen
DO 100 IY=KY1,KY2 ! added here so users can leave
DO 100 IX=KX1,KX2 ! certain points to be Neumann,
if (tmask(ix,iy,iz,iel,ifld1).eq.0) then ! if desired.
CALL NEKASGN (IX,IY,IZ,IEL)
CALL USERBC (IX,IY,IZ,IFACE,IEG)
S(IX,IY,IZ) = TEMP
endif
100 CONTINUE
RETURN
elseif (cb.eq.'o ') then
DO 101 IZ=KZ1,KZ2 ! 11/19/2010: The tmask() screen
DO 101 IY=KY1,KY2 ! added here so users can leave
DO 101 IX=KX1,KX2 ! certain points to be Neumann,
CALL NEKASGN (IX,IY,IZ,IEL)
CALL USERBC (IX,IY,IZ,IFACE,IEG)
S(IX,IY,IZ) = PA
101 CONTINUE
RETURN
ELSEIF (CB.EQ.'ms ' .OR. CB.EQ.'msi') THEN
DO 200 IZ=KZ1,KZ2
DO 200 IY=KY1,KY2
DO 200 IX=KX1,KX2
CALL NEKASGN (IX,IY,IZ,IEL)