-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathscan.f90
1475 lines (1114 loc) · 40.7 KB
/
scan.f90
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
! Part of Zinc FE package. Author: John Blackburn
module scan
! ----------------------------------------------------------------------
! Set u and ftol and variables in "common" BEFORE calling linescan, planescan
! Note that one of the args of line/planescan is type(Texpr)
! so this type is made available as well
! ----------------------------------------------------------------------
use evaluate, only : Tlist
implicit none
save
private
public u,ftol,linescan,planescan,surfint,volint,Texpr,findu,fdchk,replaceBrace
double precision, allocatable :: u(:,:,:,:)
double precision ftol
type Texpr
double precision expr
character(1000) sexpr
type(Tlist) texpr
integer iexpr
integer nbrak
integer, allocatable :: nperbrak(:) ! no. terms in each bracket
character(1000), allocatable :: sbrak(:,:) ! string for each bracket term
double precision, allocatable :: brak(:,:) ! value for each bracket term
integer, allocatable :: ibrak(:,:) ! index for each bracket term
integer, allocatable :: rbrak(:,:) ! region for each bracket term
type(Tlist), allocatable :: tbrak(:,:) ! token list for each bracket term
end type Texpr
contains
! ######################################################################
subroutine getu(i,j,k,x,y,z,ur,dur,xi,eta,mu,ifail)
! ----------------------------------------------------------------------
! For a given point r=(x,y,z) in a given element specified by its
! low indices i,j,k, calculate the value of u_ii and its derivatives
! dur(1,2)=du1/dx2 etc.
! If point not in box containing element, ifail=2
! If point is in box but not element return ifail=1
! success: ifail=0
! shares block fcnb with sub fcn
! Uses u,ftol from common block
! ----------------------------------------------------------------------
use fcnb
use shape
use geom
use common, only : nvar,imax,jmax,kmax,rnode
integer i,j,k,ifail
double precision x,y,z
double precision xi,eta,mu,ur(nvar),dur(nvar,3) ! dimensioned from common
integer, parameter :: lwa=3*(3+5)+3 ! n*(m+5)+m for iopt=1 or 2 (dnls1e)
integer iw(3),itab(8,4,3),itabsrt(8,3),icomb(2,3)
double precision rtab(8,3),wa(lwa),xvec(3),fvec(3)
double precision ul(8,nvar),dN(3) ! ul auto array
double precision w(100)
double precision jac(3,3),invjac(3,3),jacdet
integer ifirst,l,inode,jnode,knode,m,n,nprint,lw,ii,jj,info
double precision x1,y1,z1,xmin,ymin,zmin,xmax,ymax,zmax,tol
! external fcn,fcn2,fcnnag
! ----------------------------------------------------------------------
! Get element and its node u values
! put copy of point on fcnb common
! ----------------------------------------------------------------------
icomb(1,1)=i
icomb(1,2)=j
icomb(1,3)=k
icomb(2,1)=i+1
icomb(2,2)=j+1
icomb(2,3)=k+1
call getel(icomb,rnode,imax,jmax,kmax,itab,rtab)
! ----------------------------------------------------------------------
! First check if point x,y,z is inside box containing element
! If not, return will ifail=2
! ----------------------------------------------------------------------
ifirst=1
do l=1,8
x1=rtab(l,1)
y1=rtab(l,2)
z1=rtab(l,3)
if (ifirst.eq.1) then
xmax=x1
ymax=y1
zmax=z1
xmin=x1
ymin=y1
zmin=z1
ifirst=0
else
xmax=max(x1,xmax)
ymax=max(y1,ymax)
zmax=max(z1,zmax)
xmin=min(x1,xmin)
ymin=min(y1,ymin)
zmin=min(z1,zmin)
endif
enddo
if (.not.(x.le.xmax.and.x.ge.xmin.and. &
y.le.ymax.and.y.ge.ymin.and. &
z.le.zmax.and.z.ge.zmin)) then
ifail=2
return
endif
! ----------------------------------------------------------------------
! Prepare element in sorted order, get u at each local node (if allocated)
! ----------------------------------------------------------------------
call sort(itab,rtab,itabsrt,rtabsrt) ! rtabsrt passed to fcn
if (allocated(u)) then
do l=1,8
inode=itabsrt(l,1)
jnode=itabsrt(l,2)
knode=itabsrt(l,3)
do ii=1,nvar
ul(l,ii)=u(inode,jnode,knode,ii)
enddo
enddo
endif
xx=x ! passed to fcn
yy=y
zz=z
! ----------------------------------------------------------------------
! Initialise xvec, call dnls1e and recover solution xi,eta,mu
! fcn just provides fve! (iopt=1), fcn2 also jacobian (iopt=2)
! fcnnag for NAG library call
! Uncomment one of these
! Recommend use dnls1e with fcn (iopt=1) as others did not work
! correctly
! ----------------------------------------------------------------------
xvec(1)=0
xvec(2)=0
xvec(3)=0
m=3
n=3
tol=1e-6
nprint=-1
lw=100
call dnls1e(fcn, 1,m,n,xvec,fvec,tol,nprint,info,iw,wa,lwa)
! call dnls1e(fcn2,2,m,n,xvec,fvec,tol,nprint,info,iw,wa,lwa)
! call E04FYF(m,n,fcnnag,xvec,fsumsq,w,lw,iuser,user,ifail)
xi=xvec(1)
eta=xvec(2)
mu=xvec(3)
! ----------------------------------------------------------------------
! If local coords in range then call jacobian (for benefit of getdN)
! Hence prepare position and derivatives. ftol gives margin for error
! ----------------------------------------------------------------------
if ( xi .ge.-1-ftol.and.xi .le.+1+ftol.and.&
eta.ge.-1-ftol.and.eta.le.+1+ftol.and.&
mu .ge.-1-ftol.and.mu .le.+1+ftol) then
ifail=0
else
ifail=1
return
endif
call jacobian(xi,eta,mu,rtabsrt,jac,invjac,jacdet)
if (allocated(u)) then
do ii=1,nvar
ur(ii)=0
do jj=1,3
dur(ii,jj)=0
enddo
enddo
do l=1,8
call getdN(l,xi,eta,mu,invjac,dN)
do ii=1,nvar
ur(ii)=ur(ii)+ul(l,ii)*Nl(l,xi,eta,mu)
enddo
do ii=1,nvar
do jj=1,3
dur(ii,jj)=dur(ii,jj)+ul(l,ii)*dN(jj)
enddo
enddo
enddo
endif
end subroutine getu
! ######################################################################
subroutine fcn(iflag,m,n,xvec,fvec,fjac,ldfjac)
! ----------------------------------------------------------------------
! set fvec as diff between x1 for xvec=(xi,eta,mu) and xx target etc
! suitable for dnls1e with IOPT=1, ie, onlt function passed not jacobian
! ----------------------------------------------------------------------
use fcnb
use shape
double precision xvec(n),fvec(m),fjac
integer iflag,m,n,ldfjac,l
double precision xi,eta,mu,x1,y1,z1
xi=xvec(1)
eta=xvec(2)
mu=xvec(3)
x1=0
y1=0
z1=0
do l=1,8
x1=x1+rtabsrt(l,1)*Nl(l,xi,eta,mu)
y1=y1+rtabsrt(l,2)*Nl(l,xi,eta,mu)
z1=z1+rtabsrt(l,3)*Nl(l,xi,eta,mu)
enddo
fvec(1)=x1-xx
fvec(2)=y1-yy
fvec(3)=z1-zz
end subroutine fcn
! ######################################################################
subroutine fcnnag(m,n,xvec,fvec,iuser,user)
! ----------------------------------------------------------------------
! set fvec as difference between x1 for xvec=(xi,eta,mu) and xx target
! Suitable for nag function E04FYF
! ----------------------------------------------------------------------
use fcnb
use shape
integer m,n,iuser,l
double precision xvec(n),fvec(m),user
double precision xi,eta,mu,x1,y1,z1
xi=xvec(1)
eta=xvec(2)
mu=xvec(3)
x1=0
y1=0
z1=0
do l=1,8
x1=x1+rtabsrt(l,1)*Nl(l,xi,eta,mu)
y1=y1+rtabsrt(l,2)*Nl(l,xi,eta,mu)
z1=z1+rtabsrt(l,3)*Nl(l,xi,eta,mu)
enddo
fvec(1)=x1-xx
fvec(2)=y1-yy
fvec(3)=z1-zz
end subroutine fcnnag
! ######################################################################
subroutine fcn2(iflag,m,n,xvec,fvec,fjac,ldfjac)
! ----------------------------------------------------------------------
! set fvec as difference between x1 for xvec=(xi,eta,mu) and xx target
! This version also prepares jacobian at specified point
! Suitable for dnls1e IOPT=2
! calls functions Nl, dNdxi, dNdeta, dNdmu
! ----------------------------------------------------------------------
use fcnb
use shape
double precision xvec(n),fvec(m),fjac(ldfjac,n)
integer iflag,m,n,ldfjac,l,i,j
double precision xi,eta,mu,x1,y1,z1
xi=xvec(1)
eta=xvec(2)
mu=xvec(3)
if (iflag.eq.1) then
x1=0
y1=0
z1=0
do l=1,8
x1=x1+rtabsrt(l,1)*Nl(l,xi,eta,mu)
y1=y1+rtabsrt(l,2)*Nl(l,xi,eta,mu)
z1=z1+rtabsrt(l,3)*Nl(l,xi,eta,mu)
enddo
fvec(1)=x1-xx
fvec(2)=y1-yy
fvec(3)=z1-zz
else if (iflag.eq.2) then
do i=1,3
do j=1,3
fjac(i,j)=0
enddo
enddo
do i=1,3
do l=1,8
fjac(i,1)=fjac(i,1)+rtabsrt(l,i)*dNdxi(l,xi,eta,mu)
fjac(i,2)=fjac(i,2)+rtabsrt(l,i)*dNdeta(l,xi,eta,mu)
fjac(i,3)=fjac(i,3)+rtabsrt(l,i)*dNdmu(l,xi,eta,mu)
enddo
enddo
endif
end subroutine fcn2
! ######################################################################
subroutine linescan(iunit,N,x1,y1,z1,x2,y2,z2,expr)
! ----------------------------------------------------------------------
! Scan between points and write out to unit iunit. N intervals
! a little bit smart in that it tries the current element first
! Uses ijkmax,nvar from common
! ----------------------------------------------------------------------
use evaluate, only : Tlist
use common, only : nvar,imax,jmax,kmax,iregup
integer iunit,N
double precision x1,y1,z1,x2,y2,z2
type(Texpr) expr
double precision ur(nvar),dur(nvar,3) ! auto arrays
double precision xi,eta,mu,dx,dy,dz,x,y,z,pval,dist,nx,ny,nz
integer iold,jold,kold,is,irege,ios,ifail
nx=0; ny=0; nz=0
dx=(x2-x1)/N
dy=(y2-y1)/N
dz=(z2-z1)/N
iold=imax/2
jold=jmax/2
kold=kmax/2
! ----------------------------------------------------------------------
! Linescan. iold, jold, kold is last result. We will try this first
! ----------------------------------------------------------------------
do is=0,N
x=x1+is*dx
y=y1+is*dy
z=z1+is*dz
call findu(iold,jold,kold,x,y,z,ur,dur,xi,eta,mu,ifail)
! ----------------------------------------------------------------------
! If success, ifail=0 and we jumped to 1. If failure, the loop
! termination and ifail=1,2. In that case, set ur,dur,irege=0
! ----------------------------------------------------------------------
if (ifail.eq.0) then
irege=iregup(iold,jold,kold)
pval=getval(x,y,z,nx,ny,nz,ur,dur,expr,irege)
else
irege=0
pval=0
endif
dist=sqrt((x-x1)**2+(y-y1)**2+(z-z1)**2)
write (iunit,'(5e13.5,5i6)',iostat=ios) dist,x,y,z,pval,ifail,iold,jold,kold,irege
if (ios.ne.0) then
print *,'Error writing to linescan.out'
call exit(1)
endif
enddo
end subroutine linescan
! ######################################################################
subroutine planescan(iunit,inorm,pnorm,a1,a2,b1,b2,Na,Nb,expr)
! ----------------------------------------------------------------------
! Scan plane normal to inorm, 1=x, 2=y, 3=z, at pnorm
! number of intervals is Na, Nb between [a1,a2] etc
! eg, inorm=3, plane at z=pnorm, scan x=[a1,a2], y=[b1,b2]
! uses only imax etc from common
! ----------------------------------------------------------------------
use evaluate, only : Tlist
use common, only : nvar,imax,jmax,kmax,iregup
integer iunit,inorm,Na,Nb
double precision pnorm,a1,a2,b1,b2,nx,ny,nz
type(Texpr) expr
double precision ur(nvar),dur(nvar,3),r(3) ! auto arrays
double precision xi,eta,mu,x,y,z,pval
double precision da,db,a,b
integer iold,jold,kold,irege,ios
integer iac,ibc,ia,ib,ifail
nx=0; ny=0; nz=0
if (inorm.eq.1) then
iac=2
ibc=3
else if (inorm.eq.2) then
iac=1
ibc=3
else
iac=1
ibc=2
endif
da=(a2-a1)/Na
db=(b2-b1)/Nb
iold=imax/2
jold=jmax/2
kold=kmax/2
do ia=0,Na
do ib=0,Nb
a=a1+ia*da
b=b1+ib*db
r(iac)=a
r(ibc)=b
r(inorm)=pnorm
x=r(1)
y=r(2)
z=r(3)
call findu(iold,jold,kold,x,y,z,ur,dur,xi,eta,mu,ifail)
if (ifail.eq.0) then
irege=iregup(iold,jold,kold)
pval=getval(x,y,z,nx,ny,nz,ur,dur,expr,irege)
else
irege=0
pval=0
endif
write (iunit,'(4e13.5,5i6)',iostat=ios) x,y,z,pval,ifail,iold,jold,kold,irege
if (ios.ne.0) then
print *,'Error writing to planescan.out'
call exit(1)
endif
enddo
write (iunit,*)
enddo
end subroutine planescan
! ######################################################################
function surfint(ifrom,ito,expr)
! ----------------------------------------------------------------------
! Eg: surfint 2-1 = eps1*(Vx*nx+Vy*ny+Vz*nz)
! itabsrt, rtabsrt is face
! itabsrt1,2,el rtabsrt1,2,el are for element (same for irege12)
! centre12 are normally element centres, facecentre is obvious
! ul is for element, ulface is for face
! ----------------------------------------------------------------------
use common, only : imax,jmax,kmax,rnode,iregup,nvar,gauss,wt,ng_xi,ng_mu,ng_eta
use geom
use integrate
use matrices, only : lookind
use util
use shape
use iofile, only : prterr
double precision surfint
type(Texpr) expr
integer ifrom,ito
integer idir,itop,jtop,ktop,i,j,k,irege1,irege2
integer itabsrt(4,3),itabsrt1(8,3),itabsrt2(8,3),icomb(2,3),itabsrtel(8,3)
integer inode,jnode,knode,itab(8,4,3)
integer ixi,ieta,icoord,ii,jj,l,facelst(4)
double precision ul(8,nvar),ulface(4,nvar),ur(nvar),dur(nvar,3) ! auto arrays
double precision rtabsrt(4,3),rtab(8,3),rtabsrt1(8,3),rtabsrt2(8,3),rtabsrtel(8,3)
double precision facecentre(3),centre1(3),centre2(3)
double precision vec1(3),vec2(3),norm(3),norm1(3),norm2(3),norm3(3),norm4(3),test(3),eta,xi
double precision drdxi(3),drdeta(3),trailer,coord,nx,ny,nz
double precision jac(3,3),invjac(3,3),jacdet,dN(3),x,y,z
! print *,'surfint',ifrom,ito
surfint=0
do idir=1,3
if (idir.eq.1) then
itop=imax-1
jtop=jmax-1
ktop=kmax
else if (idir.eq.2) then
itop=imax-1
ktop=kmax-1
jtop=jmax
else
jtop=jmax-1
ktop=kmax-1
itop=imax
endif
! ----------------------------------------------------------------------
! Loop over all faces. get face geometry in itabsrt, rtabsrt
! get the region number and element setup in elements either side of the face
! elements stored in itabsrt12, rtabsrt12. Centres in centre12
! face centre is facecentre
! ----------------------------------------------------------------------
do i=0,itop
do j=0,jtop
do k=0,ktop
if (idir.eq.1) then
rtabsrt(1,:)=rnode(i,j,k,:)
rtabsrt(2,:)=rnode(i+1,j,k,:)
rtabsrt(3,:)=rnode(i+1,j+1,k,:)
rtabsrt(4,:)=rnode(i,j+1,k,:)
itabsrt(1,:)=[i,j,k]
itabsrt(2,:)=[i+1,j,k]
itabsrt(3,:)=[i+1,j+1,k]
itabsrt(4,:)=[i,j+1,k]
facecentre=(rtabsrt(1,:)+rtabsrt(2,:)+rtabsrt(3,:)+rtabsrt(4,:))/4
if (k.eq.kmax) then
irege1=lookind('ZMAX')
centre1=facecentre
else
irege1=iregup(i,j,k)
centre1=centre(i,j,k)
icomb(1,:)=[i,j,k]
icomb(2,:)=[i+1,j+1,k+1]
call getel(icomb,rnode,imax,jmax,kmax,itab,rtab)
call sort(itab,rtab,itabsrt1,rtabsrt1)
endif
if (k.eq.0) then
irege2=lookind('ZMIN')
centre2=facecentre
else
irege2=iregup(i,j,k-1)
centre2=centre(i,j,k-1)
icomb(1,:)=[i,j,k-1]
icomb(2,:)=[i+1,j+1,k]
call getel(icomb,rnode,imax,jmax,kmax,itab,rtab)
call sort(itab,rtab,itabsrt2,rtabsrt2)
endif
! ----------------------------------------------------------------------
else if (idir.eq.2) then
rtabsrt(1,:)=rnode(i,j,k,:)
rtabsrt(2,:)=rnode(i+1,j,k,:)
rtabsrt(3,:)=rnode(i+1,j,k+1,:)
rtabsrt(4,:)=rnode(i,j,k+1,:)
itabsrt(1,:)=[i,j,k]
itabsrt(2,:)=[i+1,j,k]
itabsrt(3,:)=[i+1,j,k+1]
itabsrt(4,:)=[i,j,k+1]
facecentre=(rtabsrt(1,:)+rtabsrt(2,:)+rtabsrt(3,:)+rtabsrt(4,:))/4
if (j.eq.jmax) then
irege1=lookind('YMAX')
centre1=facecentre
else
irege1=iregup(i,j,k)
centre1=centre(i,j,k)
icomb(1,:)=[i,j,k]
icomb(2,:)=[i+1,j+1,k+1]
call getel(icomb,rnode,imax,jmax,kmax,itab,rtab)
call sort(itab,rtab,itabsrt1,rtabsrt1)
endif
if (j.eq.0) then
irege2=lookind('YMIN')
centre2=facecentre
else
irege2=iregup(i,j-1,k)
centre2=centre(i,j-1,k)
icomb(1,:)=[i,j-1,k]
icomb(2,:)=[i+1,j,k+1]
call getel(icomb,rnode,imax,jmax,kmax,itab,rtab)
call sort(itab,rtab,itabsrt2,rtabsrt2)
endif
! ----------------------------------------------------------------------
else
rtabsrt(1,:)=rnode(i,j,k,:)
rtabsrt(2,:)=rnode(i,j+1,k,:)
rtabsrt(3,:)=rnode(i,j+1,k+1,:)
rtabsrt(4,:)=rnode(i,j,k+1,:)
itabsrt(1,:)=[i,j,k]
itabsrt(2,:)=[i,j+1,k]
itabsrt(3,:)=[i,j+1,k+1]
itabsrt(4,:)=[i,j,k+1]
facecentre=(rtabsrt(1,:)+rtabsrt(2,:)+rtabsrt(3,:)+rtabsrt(4,:))/4
if (i.eq.imax) then
irege1=lookind('XMAX')
centre1=facecentre
else
irege1=iregup(i,j,k)
centre1=centre(i,j,k)
icomb(1,:)=[i,j,k]
icomb(2,:)=[i+1,j+1,k+1]
call getel(icomb,rnode,imax,jmax,kmax,itab,rtab)
call sort(itab,rtab,itabsrt1,rtabsrt1)
endif
if (i.eq.0) then
irege2=lookind('XMIN')
centre2=facecentre
else
irege2=iregup(i-1,j,k)
centre2=centre(i-1,j,k)
icomb(1,:)=[i-1,j,k]
icomb(2,:)=[i,j+1,k+1]
call getel(icomb,rnode,imax,jmax,kmax,itab,rtab)
call sort(itab,rtab,itabsrt2,rtabsrt2)
endif
endif
! ----------------------------------------------------------------------
! If this is a face to integrate, find the unit normal vector
! discover icoord and coord indicating which face of the "to" element
! is actuall the current face. icoord = 1,2,3 => xi,eta,mu
! also get ul for relevant element
! ----------------------------------------------------------------------
if ( (irege1.eq.ifrom.and.irege2.eq.ito).or. &
(irege2.eq.ifrom.and.irege1.eq.ito)) then
! print *,'found integration face:',irege1,irege2,idir,i,j,k
vec1=rtabsrt(2,:)-rtabsrt(1,:)
vec2=rtabsrt(4,:)-rtabsrt(1,:)
norm1=cross(vec1,vec2)
vec1=rtabsrt(3,:)-rtabsrt(2,:)
vec2=rtabsrt(1,:)-rtabsrt(2,:)
norm2=cross(vec1,vec2)
vec1=rtabsrt(4,:)-rtabsrt(3,:)
vec2=rtabsrt(2,:)-rtabsrt(3,:)
norm3=cross(vec1,vec2)
vec1=rtabsrt(1,:)-rtabsrt(4,:)
vec2=rtabsrt(3,:)-rtabsrt(4,:)
norm4=cross(vec1,vec2)
if (dot_product(norm1,norm2).lt.0.or.dot_product(norm1,norm3).lt.0.or.dot_product(norm1,norm4).lt.0.or.&
dot_product(norm2,norm3).lt.0.or.dot_product(norm2,norm4).lt.0.or.dot_product(norm3,norm4).lt.0) then
call prterr('incompatible normal vectors')
endif
norm=(norm1+norm2+norm3+norm4)/4
norm=norm/mag(norm)
if (irege1.eq.ifrom.and.irege2.eq.ito) then
test=centre2-centre1
call lockface(itabsrt2,itabsrt,icoord,coord,facelst)
do l=1,8
inode=itabsrt2(l,1)
jnode=itabsrt2(l,2)
knode=itabsrt2(l,3)
do ii=1,nvar
ul(l,ii)=u(inode,jnode,knode,ii)
enddo
enddo
rtabsrtel=rtabsrt2
itabsrtel=itabsrt2
else
test=centre1-centre2
call lockface(itabsrt1,itabsrt,icoord,coord,facelst)
do l=1,8
inode=itabsrt1(l,1)
jnode=itabsrt1(l,2)
knode=itabsrt1(l,3)
do ii=1,nvar
ul(l,ii)=u(inode,jnode,knode,ii)
enddo
enddo
rtabsrtel=rtabsrt1
itabsrtel=itabsrt1
endif
if (dot_product(norm,test).lt.0) norm=-norm
! ----------------------------------------------------------------------
! Get the solution at the face nodes
! ----------------------------------------------------------------------
do l=1,4
inode=itabsrtel(facelst(l),1)
jnode=itabsrtel(facelst(l),2)
knode=itabsrtel(facelst(l),3)
do ii=1,nvar
ulface(l,ii)=u(inode,jnode,knode,ii)
enddo
enddo
! print *,'ulface',ulface,ng
! ----------------------------------------------------------------------
! Integrate
! ----------------------------------------------------------------------
do ixi=1,ng_xi
do ieta=1,ng_eta
! print *,'ixi,ieta',ixi,ieta
xi=gauss(ng_xi,ixi)
eta=gauss(ng_eta,ieta)
drdxi=0
drdeta=0
do l=1,4
drdxi(:)= drdxi(:)+ rtabsrtel(facelst(l),:)*dN2d_dxi(l,xi,eta)
drdeta(:)=drdeta(:)+rtabsrtel(facelst(l),:)*dN2d_deta(l,xi,eta)
enddo
trailer=mag(cross(drdxi,drdeta))*wt(ng_xi,ixi)*wt(ng_eta,ieta)
ur=0
do ii=1,nvar
do l=1,4
ur(ii)=ur(ii)+ulface(l,ii)*N2d(l,xi,eta)
enddo
enddo
if (icoord.eq.1) then
call jacobian(coord,xi,eta,rtabsrtel,jac,invjac,jacdet)
else if (icoord.eq.2) then
call jacobian(xi,coord,eta,rtabsrtel,jac,invjac,jacdet)
else
call jacobian(xi,eta,coord,rtabsrtel,jac,invjac,jacdet)
endif
dur(:,:)=0
do l=1,8
if (icoord.eq.1) then
call getdN(l,coord,xi,eta,invjac,dN)
else if (icoord.eq.2) then
call getdN(l,xi,coord,eta,invjac,dN)
else
call getdN(l,xi,eta,coord,invjac,dN)
endif
do ii=1,nvar
do jj=1,3
dur(ii,jj)=dur(ii,jj)+ul(l,ii)*dN(jj)
enddo
enddo
enddo
x=facecentre(1)
y=facecentre(2)
z=facecentre(3)
nx=norm(1)
ny=norm(2)
nz=norm(3)
! print *,'integrate: ',xi,eta,trailer
surfint=surfint+getval(x,y,z,nx,ny,nz,ur,dur,expr,ito)*trailer
enddo
enddo
! stop
endif
enddo
enddo
enddo
enddo
! print *,'return=',surfint
end function surfint
! ######################################################################
function volint(region,expr)
use common
use geom
use shape
use indexq, only : indQ
integer region
type(Texpr) expr
double precision volint
integer i,j,k,irege,icomb(2,3),itabsrt(8,3),itab(8,4,3)
double precision rtab(8,3),rtabsrt(8,3)
double precision ul(8,nvar),nx,ny,nz,x,y,z,ur(nvar),dur(nvar,3)
integer ixi,ieta,imu,ii,jj,inode,jnode,knode,l
double precision xi,eta,mu,trailer,dN(3),Nl1
double precision jac(3,3),invjac(3,3),jacdet
nx=0
ny=0
nz=0
volint=0
do i=0,imax-1
do j=0,jmax-1
do k=0,kmax-1
! ----------------------------------------------------------------------
! Get region number
! ----------------------------------------------------------------------
irege=iregup(i,j,k)
if (irege.eq.region) then
! ----------------------------------------------------------------------
! Get the sorted element nodes. Get local u values in ul(l,ii)
! ----------------------------------------------------------------------
icomb(1,1)=i
icomb(1,2)=j
icomb(1,3)=k
icomb(2,1)=i+1
icomb(2,2)=j+1
icomb(2,3)=k+1
call getel(icomb,rnode,imax,jmax,kmax,itab,rtab)
call sort(itab,rtab,itabsrt,rtabsrt)
x=sum(rtabsrt(:,1))/8
y=sum(rtabsrt(:,2))/8
z=sum(rtabsrt(:,3))/8
do l=1,8
inode=itabsrt(l,1)
jnode=itabsrt(l,2)
knode=itabsrt(l,3)
do ii=1,nvar
ul(l,ii)=u(inode,jnode,knode,ii)
enddo
enddo
! ----------------------------------------------------------------------
! Element integral loops
! Integrate over element, get du(ii,jj), ur(ii)
! ----------------------------------------------------------------------
do ixi=1,ng_xi
do ieta=1,ng_eta
do imu=1,ng_mu
xi=gauss(ng_xi,ixi)
eta=gauss(ng_eta,ieta)
mu=gauss(ng_mu,imu)
call jacobian(xi,eta,mu,rtabsrt,jac,invjac,jacdet)
trailer=jacdet*wt(ng_xi,ixi)*wt(ng_eta,ieta)*wt(ng_mu,imu)
dur(:,:)=0
do l=1,8
call getdN(l,xi,eta,mu,invjac,dN)
do ii=1,nvar
do jj=1,3
dur(ii,jj)=dur(ii,jj)+ul(l,ii)*dN(jj)
enddo
enddo
enddo
ur(:)=0
do l=1,8
Nl1=Nl(l,xi,eta,mu)
do ii=1,nvar
ur(ii)=ur(ii)+ul(l,ii)*Nl1
enddo
enddo
volint=volint+getval(x,y,z,nx,ny,nz,ur,dur,expr,region)*trailer
enddo
enddo
enddo
endif
enddo
enddo
enddo
end function volint
! ######################################################################
subroutine findu(iold,jold,kold,x,y,z,ur,dur,xi,eta,mu,ifail)
! ----------------------------------------------------------------------
! find ur, dur for point (x,y,z). iold,jold,kold is initial guess for element
! (bottom left corner node)
! and will be reset to the true value of the element
! ifail=0 success, ifail=1 failure
! ----------------------------------------------------------------------
use common, only : nvar,imax,jmax,kmax
integer iold,jold,kold,ifail
double precision x,y,z,ur(nvar),dur(nvar,3),xi,eta,mu
integer ijkmax,irad,i1,i2,j1,j2,k1,k2,i,j,k
logical iedge,jedge,kedge
ijkmax=max(imax,jmax,kmax)
do irad=0,ijkmax
i1=max(iold-irad,0)
i2=min(iold+irad,imax-1)