-
-
Notifications
You must be signed in to change notification settings - Fork 20
/
getflux.f
338 lines (325 loc) · 11.3 KB
/
getflux.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
!
! CalculiX - A 3-dimensional finite element program
! Copyright (C) 1998-2015 Guido Dhondt
!
! 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(version 2);
!
!
! 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., 675 Mass Ave, Cambridge, MA 02139, USA.
!
subroutine getflux(co,ntmat_,vold,
& cocon,ncocon,iset,istartset,iendset,ipkon,lakon,kon,
& ialset,ielmat,mi,flux)
!
! This is a copy of printoutface.f, modified to extract heat flux values
! at the specifice face set, to be used for coupling with preCICE
!
implicit none
!
character*8 lakonl,lakon(*)
!
integer konl(20),ifaceq(8,6),nelem,ii,i,j,i1,i2,j1,
& ncocon(2,*),k1,jj,ig,ntmat_,nope,nopes,imat,
& mint2d,ifacet(6,4),ifacew(8,5),iflag,indexe,jface,istartset(*),
& iendset(*),ipkon(*),kon(*),iset,ialset(*),
& mi(*),ielmat(mi(3),*),fidx
!
real*8 co(3,*),xl(3,20),shp(4,20),xs2(3,7),dvi,f(0:3),
& vkl(0:3,3),t(3,3),div,
& voldl(0:mi(2),20),cocon(0:6,ntmat_,*),xl2(3,8),xsj2(3),
& shp2(7,8),vold(0:mi(2),*),xi,et,xsj,temp,xi3d,et3d,ze3d,weight,
& xlocal20(3,9,6),xlocal4(3,1,4),xlocal10(3,3,4),xlocal6(3,1,5),
& xlocal15(3,4,5),xlocal8(3,4,6),xlocal8r(3,1,6),pres,
& tf(0:3),tn,tt,dd,coords(3),cond, flux(*), avgflux
!
include "gauss.f"
include "xlocal.f"
!
data ifaceq /4,3,2,1,11,10,9,12,
& 5,6,7,8,13,14,15,16,
& 1,2,6,5,9,18,13,17,
& 2,3,7,6,10,19,14,18,
& 3,4,8,7,11,20,15,19,
& 4,1,5,8,12,17,16,20/
data ifacet /1,3,2,7,6,5,
& 1,2,4,5,9,8,
& 2,3,4,6,10,9,
& 1,4,3,8,10,7/
data ifacew /1,3,2,9,8,7,0,0,
& 4,5,6,10,11,12,0,0,
& 1,2,5,4,7,14,10,13,
& 2,3,6,5,8,15,11,14,
& 4,6,3,1,12,15,9,13/
data iflag /3/
!
!
!
! initialisierung of the flux
!
f(0)=0.d0
! index for the output array
fidx=1
!
!
do jj=istartset(iset),iendset(iset)
!
jface=ialset(jj)
!
nelem=int(jface/10.d0)
ig=jface-10*nelem
lakonl=lakon(nelem)
indexe=ipkon(nelem)
imat=ielmat(1,nelem) ! what is this?
!
if(lakonl(4:4).eq.'2') then
nope=20
nopes=8
elseif(lakonl(4:4).eq.'8') then
nope=8
nopes=4
elseif(lakonl(4:5).eq.'10') then
nope=10
nopes=6
elseif(lakonl(4:4).eq.'4') then
nope=4
nopes=3
elseif(lakonl(4:5).eq.'15') then
nope=15
elseif(lakonl(4:4).eq.'6') then
nope=6
endif
!
if(lakonl(4:5).eq.'8R') then
mint2d=1
elseif((lakonl(4:4).eq.'8').or.(lakonl(4:6).eq.'20R'))
& then
if((lakonl(7:7).eq.'A').or.(lakonl(7:7).eq.'S').or.
& (lakonl(7:7).eq.'E')) then
mint2d=2
else
mint2d=4
endif
elseif(lakonl(4:4).eq.'2') then
mint2d=9
elseif(lakonl(4:5).eq.'10') then
mint2d=3
elseif(lakonl(4:4).eq.'4') then
mint2d=1
endif
!
! local topology
!
do i=1,nope
konl(i)=kon(indexe+i)
enddo
!
! computation of the coordinates of the local nodes
!
do i=1,nope
do j=1,3
xl(j,i)=co(j,konl(i))
enddo
enddo
!
! temperature, velocity and auxiliary variables
! (rho*energy density, rho*velocity and rho)
!
do i1=1,nope
do i2=0,mi(2)
voldl(i2,i1)=vold(i2,konl(i1))
enddo
enddo
!
! treatment of wedge faces
!
if(lakonl(4:4).eq.'6') then
mint2d=1
if(ig.le.2) then
nopes=3
else
nopes=4
endif
endif
if(lakonl(4:5).eq.'15') then
if(ig.le.2) then
mint2d=3
nopes=6
else
mint2d=4
nopes=8
endif
endif
!
if((nope.eq.20).or.(nope.eq.8)) then
do i=1,nopes
do j=1,3
xl2(j,i)=co(j,konl(ifaceq(i,ig)))
enddo
enddo
elseif((nope.eq.10).or.(nope.eq.4)) then
do i=1,nopes
do j=1,3
xl2(j,i)=co(j,konl(ifacet(i,ig)))
enddo
enddo
else
do i=1,nopes
do j=1,3
xl2(j,i)=co(j,konl(ifacew(i,ig)))
enddo
enddo
endif
!
! flux for each face will be the average flux
! of all integration points
avgflux=0
do i=1,mint2d
!
! local coordinates of the surface integration
! point within the surface local coordinate system
!
if((lakonl(4:5).eq.'8R').or.
& ((lakonl(4:4).eq.'6').and.(nopes.eq.4))) then
xi=gauss2d1(1,i)
et=gauss2d1(2,i)
weight=weight2d1(i)
elseif((lakonl(4:4).eq.'8').or.
& (lakonl(4:6).eq.'20R').or.
& ((lakonl(4:5).eq.'15').and.(nopes.eq.8))) then
xi=gauss2d2(1,i)
et=gauss2d2(2,i)
weight=weight2d2(i)
elseif(lakonl(4:4).eq.'2') then
xi=gauss2d3(1,i)
et=gauss2d3(2,i)
weight=weight2d3(i)
elseif((lakonl(4:5).eq.'10').or.
& ((lakonl(4:5).eq.'15').and.(nopes.eq.6))) then
xi=gauss2d5(1,i)
et=gauss2d5(2,i)
weight=weight2d5(i)
elseif((lakonl(4:4).eq.'4').or.
& ((lakonl(4:4).eq.'6').and.(nopes.eq.3))) then
xi=gauss2d4(1,i)
et=gauss2d4(2,i)
weight=weight2d4(i)
endif
!
! local surface normal
!
if(nopes.eq.8) then
call shape8q(xi,et,xl2,xsj2,xs2,shp2,iflag)
elseif(nopes.eq.4) then
call shape4q(xi,et,xl2,xsj2,xs2,shp2,iflag)
elseif(nopes.eq.6) then
call shape6tri(xi,et,xl2,xsj2,xs2,shp2,iflag)
else
call shape3tri(xi,et,xl2,xsj2,xs2,shp2,iflag)
endif
!
! global coordinates of the integration point
!
do j1=1,3
coords(j1)=0.d0
do i1=1,nopes
coords(j1)=coords(j1)+shp2(4,i1)*xl2(j1,i1)
enddo
enddo
!
! local coordinates of the surface integration
! point within the element local coordinate system
!
if(lakonl(4:5).eq.'8R') then
xi3d=xlocal8r(1,i,ig)
et3d=xlocal8r(2,i,ig)
ze3d=xlocal8r(3,i,ig)
call shape8h(xi3d,et3d,ze3d,xl,xsj,shp,iflag)
elseif(lakonl(4:4).eq.'8') then
xi3d=xlocal8(1,i,ig)
et3d=xlocal8(2,i,ig)
ze3d=xlocal8(3,i,ig)
call shape8h(xi3d,et3d,ze3d,xl,xsj,shp,iflag)
elseif(lakonl(4:6).eq.'20R') then
xi3d=xlocal8(1,i,ig)
et3d=xlocal8(2,i,ig)
ze3d=xlocal8(3,i,ig)
call shape20h(xi3d,et3d,ze3d,xl,xsj,shp,iflag)
elseif(lakonl(4:4).eq.'2') then
xi3d=xlocal20(1,i,ig)
et3d=xlocal20(2,i,ig)
ze3d=xlocal20(3,i,ig)
call shape20h(xi3d,et3d,ze3d,xl,xsj,shp,iflag)
elseif(lakonl(4:5).eq.'10') then
xi3d=xlocal10(1,i,ig)
et3d=xlocal10(2,i,ig)
ze3d=xlocal10(3,i,ig)
call shape10tet(xi3d,et3d,ze3d,xl,xsj,shp,iflag)
elseif(lakonl(4:4).eq.'4') then
xi3d=xlocal4(1,i,ig)
et3d=xlocal4(2,i,ig)
ze3d=xlocal4(3,i,ig)
call shape4tet(xi3d,et3d,ze3d,xl,xsj,shp,iflag)
elseif(lakonl(4:5).eq.'15') then
xi3d=xlocal15(1,i,ig)
et3d=xlocal15(2,i,ig)
ze3d=xlocal15(3,i,ig)
call shape15w(xi3d,et3d,ze3d,xl,xsj,shp,iflag)
elseif(lakonl(4:4).eq.'6') then
xi3d=xlocal6(1,i,ig)
et3d=xlocal6(2,i,ig)
ze3d=xlocal6(3,i,ig)
call shape6w(xi3d,et3d,ze3d,xl,xsj,shp,iflag)
endif
!
! calculating of
! the temperature temp
! in the integration point
!
temp=0.d0
do j1=1,3
vkl(0,j1)=0.d0
enddo
do i1=1,nope
temp=temp+shp(4,i1)*voldl(0,i1)
do k1=1,3
vkl(0,k1)=vkl(0,k1)+shp(k1,i1)*voldl(0,i1)
enddo
enddo
!
! material data (conductivity)
!
call materialdata_cond(imat,ntmat_,temp,cocon,
& ncocon,cond)
!
! determining the stress
!
dd=dsqrt(xsj2(1)*xsj2(1)+xsj2(2)*xsj2(2)+
& xsj2(3)*xsj2(3))
!
tf(0)=-cond*(vkl(0,1)*xsj2(1)+
& vkl(0,2)*xsj2(2)+
& vkl(0,3)*xsj2(3))
f(0)=f(0)+tf(0)*weight
tf(0)=tf(0)/dd
! print *, cond, temp, vkl(0,2), xsj2(2), tf(0), dd
! print *, temp, tf(0), dd, temp-tf(0)*dd/cond
avgflux=avgflux+tf(0)
!
enddo
flux(fidx)=avgflux/mint2d
! print *, fidx, flux(fidx)
fidx=fidx+1
enddo
!
!
return
end