forked from Nek5000/NekExamples
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathhemi.usr
580 lines (500 loc) · 16 KB
/
hemi.usr
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
c-----------------------------------------------------------------------
c
c This example illustrates particle tracking being used in a "bubble wire"
c mode in which particles are periodicly released along a horizontal
c line segment upstream of a hemispherical roughness element [1,2].
c
c PARTICLE TRACKING
c
c Currently, this is set up to support only a few bubble wires,
c in order to illustrate basic particle tracking.
c
c CONFIGURATION:
c
c The hemisphere of radius R=0.5 is placed onto a thin cylinder of
c radius R and height dz = 0.1R, which is placed onto a smooth
c plate. The depth, width and length of the channel are H=6.5R,
c B=12.8R and L=36.4R. The boundary conditions (BC) are a Dirichlet
c BC as a prescribed Blasius profile at the inlet, outflow BC
c at the outlet and stress free symmetry conditions at the lateral
c boundaries and at the channel surface.
c
c The Blasius profile is approximated with the following formula:
c
c u = U0 sin ( pi y / 2d )
c
c d=1.2R is the boundary layer thickness
c U0 =1.0 is the free stream velocity.
c
c Suggested Reynolds numbers (Re=UR/nu): Re=100, Re=450, Re=700
c
c Particles can be plotted with the attached matlab script, "plotit.m"
c
c A visit script will be provided soon.
c
c
c
c [1] M.S. Acalar and C.R. Smith, "A study of hairpin vortices in a laminar
c boundary layer: Part 1, hairpin vortices generated by a hemisphere
c protuberance," J. Fluid Mech. 175, pp. 1-41, 1987
c
c [2] H.M. Tufo, P.F. Fischer, M.E. Papka and K. Blom, "Numerical simulation
c and immersive visualization of hairpin vortices," Proc. of the ACM/IEEE
c SC99 Conf. on High Performance Networking and Computing, IEEE Computer Soc.,
c 1999.
c
c
c-----------------------------------------------------------------------
subroutine uservp (ix,iy,iz,eg)
include 'SIZE'
include 'TOTAL'
include 'NEKUSE'
real visc_inf
save visc_inf
data visc_inf /0/
integer e,f,eg
e = gllel(eg)
udiff =0
utrans=0
return
end
c-----------------------------------------------------------------------
subroutine set_prof
include 'SIZE'
include 'TOTAL'
common /uprof/ up(lx1,ly1,lz1,lelv)
integer icalld
save icalld
data icalld /0/
if (icalld.ne.0) return
icalld = icalld + 1
c Compute the Blasius profile
u0 = 1.0
visc =param(2)
c delta=param(71)
c if (delta.eq.0.0) DELTA= 0.58
rad = 0.5
delta = 1.2*rad ! Thorsten
re = 1./visc
hk = 0.55
uk = blasius_sin(hk,delta,u0)
rk = uk*hk/visc
write(6,1) nid,delta,re,rk
1 format(i6,' Delta,Re,Rek:',1p3e15.4)
ntot=nx1*ny1*nz1*nelv
do i=1,ntot
up(i,1,1,1)=blasius_sin(zm1(i,1,1,1),delta,u0)
enddo
return
end
c-----------------------------------------------------------------------
function blasius(y,delta,u)
c Return the velocity at a given y value for specified delta and U.
integer icalld
save icalld
data icalld /0/
dimension u0(45),y0(45),work(45)
save u0 ,y0 ,work
data u0 / 0.00000 , 0.06641 , 0.13277 , 0.19894 , 0.26471
$ , 0.32979 , 0.39378 , 0.45627 , 0.51676 , 0.57477
$ , 0.62977 , 0.68132 , 0.72899 , 0.77246 , 0.81152
$ , 0.84605 , 0.87609 , 0.90177 , 0.92333 , 0.94112
$ , 0.95552 , 0.96696 , 0.97587 , 0.98269 , 0.98779
$ , 0.99155 , 0.99425 , 0.99616 , 0.99748 , 0.99838
$ , 0.99898 , 0.99937 , 0.99961 , 0.99977 , 0.99987
$ , 0.99992 , 0.99996 , 0.99998 , 0.99999 , 1.00000
$ , 1.00000 , 1.00000 , 1.00000 , 1.00000 , 1.00000 /
if (icalld.eq.0) then
c Initialize Blasius profile and spline fitting routine.
icalld=1
do 10 i=1,45
y0(i)=float(i-1)/5.0
10 continue
call spline(y0,u0,45,work)
endif
eta=5.0*y/delta
if (eta.gt.8.5) then
blasius=u
else
call splint(y0,u0,work,45,eta,vel)
blasius=vel*u
endif
return
end
c-----------------------------------------------------------------------
function blasius_sin(y,delta,U)
c Return the velocity at a given y value for specified delta and U.
real pi,one
save pi,one
data pi,one /0.,1./
if (pi.eq.0) pi=4*atan(one)
arg = .5*pi*y/delta
blasius_sin = sin(arg)
if (arg.gt.0.5*pi) blasius_sin = 1.
blasius_sin = U*blasius_sin
return
end
c-----------------------------------------------------------------------
subroutine spline(x,y,n,y2)
parameter (nmax=100)
dimension x(n),y(n),y2(n),u(nmax)
if (n.gt.nmax) then
write(6,11) n,nmax
11 FORMAT(2X,'ERROR: Attempt to fit a spline with',I5
$ ,'greater than',I4,' points.'
$ ,/,2X,'Recompile routine SPLINE.')
call exitti('routine spline fail$',n)
endif
y2(1)=0.0
u(1) =0.0
do 10 i=2,n-1
ir=i+1
il=i-1
sig=(x(i)-x(il))/(x(ir)-x(il))
p=sig*y2(il)+2.
y2(i)=(sig-1.)/p
u(i)= ( 6.*
$ ( (y(ir)-y(i))/(x(ir)-x(i))-(y(i)-y(il))/ (x(i)-x(il) ) )
$ / (x(ir)-x(il))
$ - sig*u(il) )/p
10 continue
qn=0.0
un=0.0
y2(n)=(un-qn*u(n-1))/(qn*y2(n-1)+1.)
do 20 k=n-1,1,-1
y2(k)=y2(k)*y2(k+1)+u(k)
20 continue
return
end
c-----------------------------------------------------------------------
subroutine splint(xa,ya,y2a,n,x,y)
c p. 88-89, numerical recipes
real xa(n),ya(n),y2a(n)
klo=1
khi=n
1 if ((khi-klo).gt.1) then
k=(khi+klo)/2
if (xa(k).gt.x) then
khi=k
else
klo=k
endif
goto 1
endif
h=xa(khi)-xa(klo)
if (h.eq.0) then
write(6,*) xa(khi), 'splint failure',khi
return
endif
a=(xa(khi)-x)/h
b=(x-xa(klo))/h
y=a*ya(klo)+b*ya(khi)+
$ ((a**3-a)*y2a(klo)+(b**3-b)*y2a(khi))*(h**2)/6.
return
end
c-----------------------------------------------------------------------
subroutine userf (ix,iy,iz,iel)
include 'SIZE'
include 'TSTEP'
c include 'TOTAL'
include 'NEKUSE'
ffx = 0.0
ffy = 0.0
ffz = 0.0
return
end
c-----------------------------------------------------------------------
subroutine userq (ix,iy,iz,iel)
include 'SIZE'
include 'TOTAL'
include 'NEKUSE'
qvol = 0.0
source = 0.0
return
end
c-----------------------------------------------------------------------
subroutine userchk
include 'SIZE'
include 'TOTAL'
common /rparts/ pts(ldim,lpart),vel(ldim,2:3,lpart)
common /iparts/ npart,partid(lpart)
integer partid
ifto = .true.
if (istep.eq.0.or.mod(istep,iostep).eq.0)
$ call lambda2(t) ! Put lambda2 into temperature field
c NOTE: The following outpost call coordinates the .f0000 numbering
c with the particle file numbering.
c
c This coordination is important when synching frames for a movie.
ipstep = iostep ! Particle output coordinated with iostep
ipstep = iostep/4
ifxyo = .true.
if (istep.gt.iostep) ifxyo = .false.
if (istep.eq.0) call outpost(vx,vy,vz,pr,t,' ') ! Part. coordination
call my_particle_generator(ipstep) ! Particle injection
n=nx1*ny1*nz1*nelt
umax = glamax(vx,n)
vmax = glamax(vy,n)
wmin = glmin(vz,n)
wmax = glmax(vz,n)
wmax = glmax(vz,n)
if (nid.eq.0) write(6,1) istep,time,umax,vmax,wmin,wmax
1 format(i9,1p5e12.4,' wmax')
return
end
c-----------------------------------------------------------------------
subroutine userbc (ix,iy,iz,iside,ieg)
include 'SIZE'
include 'TOTAL'
include 'NEKUSE'
common /uprof/ up(lx1,ly1,lz1,lelv)
integer icalld
save icalld
data icalld /0/
if (icalld.eq.0) call set_prof
icalld = icalld + 1
iel = gllel(ieg)
ux=up(ix,iy,iz,iel)
uy=0.0
uz=0.0
temp=1.0
return
end
c-----------------------------------------------------------------------
subroutine useric (ix,iy,iz,ieg)
include 'SIZE'
include 'TOTAL'
include 'NEKUSE'
common /uprof/ up(lx1,ly1,lz1,lelv)
integer icalld
save icalld
data icalld /0/
if (icalld.eq.0) call set_prof
icalld = icalld + 1
iel = gllel(ieg)
UX=up(ix,iy,iz,iel)
UY=0.0
UZ=0.0
TEMP=1.0
return
end
c-----------------------------------------------------------------------
subroutine usrdat3
return
end
c-----------------------------------------------------------------------
subroutine usrdat2
include 'SIZE'
include 'TOTAL'
c param(66) = 4
c param(67) = 4
return
end
c-----------------------------------------------------------------------
subroutine usrdat
return
end
c-----------------------------------------------------------------------
subroutine interp_v(uvw,xyz,n)
c
c evaluate velocity for list of points xyz
c
c Note: -- modify
c intpts to get rid off ' WARNING: point on boundary or ...'
include 'SIZE'
include 'TOTAL'
real uvw(ldim,n),xyz(ldim,n)
logical ifjac,ifpts
parameter(nmax=lpart,nfldmax=ldim)
common /rv_intp/ pts(ldim*nmax)
common /iv_intp/ ihandle
common /outtmp/ wrk(lx1*ly1*lz1*lelt,nfldmax)
integer icalld,e
save icalld
data icalld /0/
nxyz = nx1*ny1*nz1
ntot = nxyz*nelt
if (n.gt.nmax) call exitti ('ABORT: interp_v() n > nmax!$',n)
if (nelgt.ne.nelgv) call exitti
$ ('ABORT: interp_v() nelgt.ne.nelgv not yet supported!$',nelgv)
do i=1,n ! ? not moving -> save?
pts(i) = xyz(1,i)
pts(i + n) = xyz(2,i)
if (if3d) pts(i + n*2) = xyz(3,i)
enddo
if (icalld.eq.0) then ! interpolation setup !? intpts_done(ih_intp_v)?
icalld = 1
tolin = 1.e-8
call intpts_setup(tolin,ihandle)
endif
nflds = ndim ! number of fields to interpolate
! pack working array
call opcopy(wrk(1,1),wrk(1,2),wrk(1,3),vx,vy,vz)
! interpolate
ifjac = .true. ! output transpose (of Jacobian)
ifpts = .true. ! find points
call intpts(wrk,nflds,pts,n,uvw,ifjac,ifpts,ihandle) ! copy array instead?
return
end
c-----------------------------------------------------------------------
subroutine particle_advect_std(x,vh,partv,npart)
c Lagrangian particle advection
include 'SIZE'
include 'TOTAL'
real x(ldim,lpart),vh(ldim,2:3,lpart),partv(lpart)
common /scruz/ u1(ldim,lpart)
common /padvc/ xmx(3,0:2)
if (istep.eq.0) then ! AB1
call rzero(vh,3*ndim*npart)
c1 = 1.
c2 = 0.
c3 = 0.
elseif (istep.eq.1) then ! AB2
c1 = 3
c2 = -1.
c3 = 0
c1 = c1/2.
c2 = c2/2.
else ! AB3
c1 = 23.
c2 = -16.
c3 = 5.
c1 = c1/12.
c2 = c2/12.
c3 = c3/12
endif
call interp_v(u1,x,npart)
do i=1,npart
do k=1,ndim
C Update particle position and history
x(k,i) = x(k,i)
$ + dt*(c1*u1(k,i) + c2*vh(k,2,i) + c3*vh(k,3,i))
C Update particle and fluid velocity history
vh(k,3,i) = vh(k,2,i)
vh(k,2,i) = u1(k,i)
enddo
enddo
return
end
c----------------------------------------------------------------------
subroutine particle_out (x,partid,partv,npart,ipstep)
include 'SIZE'
include 'TOTAL'
real x(ldim,lpart),partv(lpart)
integer partid(lpart)
common /scrns/ x_tmp(ldim+1,lpart),work(ldim+1,lpart)
character*128 fname
integer icalld
save icalld
data icalld /0/
if (mod(istep,ipstep).ne.0) return
icalld = icalld+1
if (nid.eq.0) then
write(fname,1) icalld
1 format('part',i5.5,'.3D')
open(unit=72,file=fname)
endif
min_points = iglmin(partid,npart)
max_points = iglmax(partid,npart)
n_active = max_points - min_points + 1
npass = n_active / lpart
if (n_active.gt.npass*lpart) npass = npass+1
ilast=min_points-1
do ipass=1,npass
mpart = min(lpart,max_points-ilast)
i0 = ilast
i1 = i0 + mpart
ilast = i1
call rzero(x_tmp,(ldim+1)*lpart)
do ii=1,npart
if (i0.lt.partid(ii).and.partid(ii).le.i1) then
i = partid(ii)-i0
call copy(x_tmp(1,i),x(1,ii),ldim) ! Coordinates
x_tmp(ldim+1,i) = partv(ii) ! Store value here
endif
enddo
call gop(x_tmp,work,'+ ',(ldim+1)*lpart)
if (nid.eq.0) write(72,2)((x_tmp(k,i),k=1,ldim+1),i=1,mpart)
2 format(1p4e17.9)
enddo
if (nid.eq.0) close(72)
return
end
C-----------------------------------------------------------------------
subroutine my_particle_generator(ipstep) ! Particle injection
include 'SIZE'
include 'TOTAL'
include 'mpif.h'
common /rparts/ pts(ldim,lpart),vel(ldim,2:3,lpart),partv(lpart)
common /iparts/ npart,partid(lpart)
integer partid
real ptime ! Time the particle motion
save ptime
data ptime /0./
ptime0 = dnekclock_sync()
call particle_init (pts,partid,partv,npart,ipstep)
call particle_advect_std (pts,vel,partv,npart)
call particle_out (pts,partid,partv,npart,ipstep)
ptime1 = dnekclock_sync() ! Track computational cost of particles
dptime = (ptime1-ptime0)
ptime = ptime + dptime
if (mod(istep,ipstep).eq.0) then
npmin = iglmin(npart,1)
npmax = iglmax(npart,1)
nptot = iglsum(npart,1)
npmav = nptot/np
ptpp = dptime / nptot ! particle-time per point
if (nid.eq.0) write(6,1)
$ istep,npmin,npmav,npmax,nptot,ptpp,dptime,ptime,time
1 format(4i7,i9,1p4e15.7,' ptime')
endif
return
end
c-----------------------------------------------------------------------
subroutine particle_init(x,partid,partv,npart,ipstep)
c
c This version does continuous injection and overwrites old entries
c
include 'SIZE'
include 'TOTAL'
real x(ldim,lpart),partv(lpart)
integer partid(lpart)
integer lcount,icount
save lcount,icount
data lcount,icount /0,0/
if (mod(istep,ipstep).ne.0) return
llpart = lpart
dx = 0.05
dy = .025
ds = (1-dy)/100
k0 = 0
k = icount ! icount = total # particles emitted
l = lcount ! Index into local particle array
nw = 100
do iline = 1,4 ! 4 lines at different heights
xp = -1.0
zp = .05 + (iline-1)*0.075
do ipart = 1,nw ! nw points on a wire
dy = 0.6*2/(nw-1)
yp = -0.6 + dy*(ipart-1)
if (mod(k,np).eq.nid) then ! Inject particle _only_ on this processor
l=l+1 ! local count
if (l.gt.llpart) l=1 ! Retire old particles by overwriting their data
x(1,l) = xp
x(2,l) = yp
x(3,l) = zp
partid(l) = k+1
c partv(l) = time ! A simple coloring scheme
c partv(l) = k0 ! Another simple coloring scheme
partv(l) = iline ! color by line
endif
k = k+1 ! Total count
enddo
k0= k0+1 ! Count for this release (for coloring only)
enddo
icount = k
lcount = l
npart = max(npart,lcount)
return
end
c-----------------------------------------------------------------------