source: CIVL/examples/fortran/nek5000/core/cmt/artvisc.f

main
Last change on this file was ea777aa, checked in by Alex Wilton <awilton@…>, 3 years ago

Moved examples, include, build_default.properties, common.xml, and README out from dev.civl.com into the root of the repo.

git-svn-id: svn://vsl.cis.udel.edu/civl/trunk@5704 fb995dde-84ed-4084-dfe6-e5aef3e2452c

  • Property mode set to 100755
File size: 17.1 KB
Line 
1 subroutine compute_entropy(s)
2! computes entropy at istep and pushes the stack down for previous
3! steps needed to compute ds/dt via finite difference (for now).
4! hardcoded for Burgers equation. More later when folded into CMT-nek
5! for Burgers, s=energy=1/2 U^2
6 include 'SIZE'
7 include 'TOTAL' ! tlag is lurking. be careful
8 include 'CMTDATA'
9! I've always seen lorder=3, but I still need three steps
10! s(:, 1, 1) ! entropy at current step
11! s(:, 2, 1) ! entropy at step n-1
12! s(:, 1, 2) ! entropy at step n-2
13 real s(lx1*ly1*lz1*lelt,lorder-1,*)
14 real ntol
15 integer e
16 data icalld /0/
17 save icalld
18
19 n=lx1*ly1*lz1
20 ntot=n*nelt
21 ntol=1.0e-10
22
23 if (icalld .eq. 0) then
24 if (nio .eq. 0) write(6,*) 'zeroing out entropy stack',istep
25 icalld=1
26 call rzero(s,ntot)
27 call rzero(s(1,1,2),ntot) ! s_{n-1}
28 call rzero(s(1,2,1),ntot) ! s_n
29 endif
30
31! compute the current entropy. This actually needs to go back in the
32! usr file because it's EOS-dependent
33 rgam=rgasref/(gmaref-1.0)
34 do i=1,ntot
35 rho=max(vtrans(i,1,1,1,irho),ntol)
36 s(i,1,1)=rgam*rho*log(pr(i,1,1,1)/(rho**gmaref))
37 enddo
38
39 if (stage .eq. 1) then
40! push the stack
41 call copy(s(1,1,2),s(1,2,1),ntot) ! s_{n-1}=s_n
42 call copy(s(1,2,1),s(1,1,1),ntot) ! fill s_n
43 endif
44
45 return
46 end
47
48!-----------------------------------------------------------------------
49
50 subroutine entropy_viscosity
51 include 'SIZE'
52 include 'TOTAL'
53 include 'CMTDATA'
54 parameter (lxyz=lx1*ly1*lz1)
55 common /scrns/ scrent(lxyz,lelt)
56 integer e
57 character*132 deathmessage
58
59 pi=4.0*atan(1.0)
60
61 n=lx1*ly1*lz1
62 ntot=n*nelt
63
64! entropy at this and lorder prior steps
65 call compute_entropy(tlag)
66! compute maxval(|S-<S>|)
67 savg = glsc2(tlag,bm1,ntot)
68 savg = -savg/volvm1
69 call cadd2(scrent,tlag,savg,ntot)
70 maxdiff = glamax(scrent,ntot)
71 if (maxdiff.le.0.0) then
72 write(6,*) 'zero maxdiff usually means NAN$'
73! write(deathmessage,*) 'zero maxdiff usually means NAN$'
74! call exittr(deathmessage,maxdiff,istep)
75 endif
76 call entropy_residual(tlag) ! fill res2
77 call copy(res2(1,1,1,1,2),res2,ntot) ! raw residual in res2
78 call wavevisc(t(1,1,1,1,3))
79 call resvisc(res2) ! overwrite res2
80 call evmsmooth(res2,t(1,1,1,1,3),.true.) ! endpoints=.false.
81 ! is intended to
82 ! preserve face states,
83 ! but this is easier to
84 ! test 1D
85! call evmsmooth(res2,t(1,1,1,1,3),.true.) ! And again.
86 call dsavg(res2) ! you DEFINITELY don't want a min here
87
88
89 return
90 end
91
92!-----------------------------------------------------------------------
93
94 subroutine piecewiseAV(shock_detector)
95 include 'SIZE'
96 include 'TOTAL'
97 include 'CMTDATA'
98 integer e
99 character*132 deathmessage
100 common /scrvh/ avmask(lelt)
101 real avmask
102 external shock_detector
103
104 pi=4.0*atan(1.0)
105
106 nxyz=lx1*ly1*lz1
107 ntot=nxyz*nelt
108
109! toggle shock detector with AV application Lv, See & Ihme (2016) JCP 322
110
111 if (time4av) then
112 call shock_detector(avmask)
113
114! old nu_max
115 call wavevisc(t(1,1,1,1,3))
116! diagnostic
117! if (stage.gt.0)then
118! do e=1,nelt
119!! do i=1,nxyz
120! write(stage*100+nid,*)xm1(i,1,1,e),ym1(i,1,1,e),t(i,1,1,e,3)
121! enddo
122! enddo
123! endif
124 do e=1,nelt
125 call cmult(t(1,1,1,e,3),avmask(e),nxyz)
126 if (avmask(e).ne.1.0) write(6,*) 'duh sir'
127 enddo
128
129 else
130 call rzero(t(1,1,1,1,3),ntot)
131! call shock_detector(avmask)
132 endif
133
134! diagnostic
135! if (stage.gt.0)then
136! do e=1,nelt
137! do i=1,nxyz
138! write(stage*1000+nid,'(5e15.7)') xm1(i,1,1,e),ym1(i,1,1,e),
139!! > t(i,1,1,e,3),avmask(e),epsebdg(e)
140! enddo
141! enddo
142! endif
143
144! call max_to_trilin(t(1,1,1,1,3))
145
146 return
147 end
148
149c-----------------------------------------------------------------------
150
151 subroutine entropy_residual(s)
152! COMPUTE R=ds/dt + div F(s) for entropy pair s and (hardcoded) F
153! and store its norm functional thing |R| in res2 where it will
154! provide artificial viscosity according to the code in
155! entropy_viscosity and the method of Guermond
156 include 'SIZE'
157 include 'TOTAL' ! tlag lurks
158 include 'CMTDATA'
159 integer e
160 real dsdtcoef(3) ! put this in /TIMESTEPCOEF/ someday
161 data dsdtcoef /1.0,1.0,0.5/
162 real s(lx1*ly1*lz1*lelt,lorder-1,ldimt) ! because it's really tlag
163
164 n=lx1*ly1*lz1
165 ntot=n*nelt
166
167 if (istep .eq. 1) return
168 rdt=1.0/(dsdtcoef(stage)*dt)
169 if (stage .eq. 1) then ! THE MOST SABOTAGEABLE PART OF THE CODE
170 call sub3(res2,s(1,1,1),s(1,1,2),ntot) ! EVALUATE s_n-s_{n-1}
171 else
172 call sub3(res2,s(1,1,1),s(1,2,1),ntot) ! EVALUATE s^{(stage)}-s_n
173 endif
174 call cmult(res2,rdt,ntot)
175! res2=ds/dt. now,
176
177!-----------------------------------------------------------------------
178! cons approach: strong-conservation form of flux divergence in entropy residual
179!-----------------------------------------------------------------------
180! get around to expanding totalh to store fluxes for whole fields and
181! properly vectorize evaluate_*_h and flux_div_integral
182 do e=1,nelt
183 call evaluate_entropy_flux(e) ! diffh. zero it before diffusion
184 call flux_div_mini(e) ! into res2 it goes.
185 enddo
186
187 return
188 end
189
190c-----------------------------------------------------------------------
191
192 subroutine evaluate_entropy_flux(e)
193! entropy flux function for entropy residual.
194! just vel*s for now
195 include 'SIZE'
196 include 'SOLN'
197 include 'INPUT'
198 include 'CMTDATA'
199 integer e
200
201 call rzero(totalh,3*lxd*lyd*lzd)
202 n=lx1*ly1*lz1
203
204 call col3(totalh(1,1),vx(1,1,1,e),tlag(1,1,1,e,1,1),n)
205 call col3(totalh(1,2),vy(1,1,1,e),tlag(1,1,1,e,1,1),n)
206 if (if3d) call col3(totalh(1,3),vz(1,1,1,e),tlag(1,1,1,e,1,1),n)
207
208 return
209 end
210
211c-----------------------------------------------------------------------
212
213 subroutine flux_div_mini(e)
214 include 'SIZE'
215 include 'INPUT'
216 include 'GEOM'
217 include 'DXYZ'
218 include 'SOLN'
219 include 'CMTDATA'
220 parameter (ldd=lx1*ly1*lz1)
221 parameter (ldg=lx1**3,lwkd=2*ldg)
222 common /ctmp1/ ur(ldd),us(ldd),ut(ldd),ju(ldd),ud(ldd),tu(ldd)
223 real ju
224
225 integer e
226
227 nrstd=ldd
228 nxyz=lx1*ly1*lz1
229 mdm1=lx1-1
230
231 call rzero(ur,nrstd)
232 call rzero(us,nrstd)
233 call rzero(ut,nrstd)
234 call rzero(ud,nrstd)
235 call rzero(tu,nrstd)
236
237 if (if3d) then
238 call local_grad3(ur,us,ut,totalh(1,1),mdm1,1,dxm1,dxtm1)
239 do i=1,nxyz
240 ud(i) = jacmi(i,e)*(rxm1(i,1,1,e)*ur(i)
241 $ + sxm1(i,1,1,e)*us(i)
242 $ + txm1(i,1,1,e)*ut(i))
243 enddo
244 call local_grad3(ur,us,ut,totalh(1,2),mdm1,1,dxm1,dxtm1)
245 do i=1,nxyz ! confirmed to have no effect in 1D
246 ud(i)=ud(i)+jacmi(i,e)*(rym1(i,1,1,e)*ur(i)
247 $ + sym1(i,1,1,e)*us(i)
248 $ + tym1(i,1,1,e)*ut(i))
249 enddo
250 call local_grad3(ur,us,ut,totalh(1,3),mdm1,1,dxm1,dxtm1)
251 do i=1,nxyz ! confirmed to have no effect in 1D
252 ud(i)=ud(i)+jacmi(i,e)*(rzm1(i,1,1,e)*ur(i)
253 $ + szm1(i,1,1,e)*us(i)
254 $ + tzm1(i,1,1,e)*ut(i))
255 enddo
256 else
257 call local_grad2(ur,us,totalh(1,1),mdm1,1,dxm1,dxtm1)
258 do i=1,nxyz
259 ud(i) = jacmi(i,e)*(rxm1(i,1,1,e)*ur(i)
260 $ + sxm1(i,1,1,e)*us(i))
261 enddo
262 call local_grad2(ur,us,totalh(1,2),mdm1,1,dxm1,dxtm1)
263 do i=1,nxyz
264 ud(i)=ud(i)+jacmi(i,e)*(rym1(i,1,1,e)*ur(i)
265 $ + sym1(i,1,1,e)*us(i))
266 enddo
267 endif
268 call add2(res2(1,1,1,e,1),ud,nxyz)
269
270 return
271 end
272
273!-----------------------------------------------------------------------
274
275 subroutine resvisc(residual)
276! smooth residual-based entropy visc, defined by Guermond, Popov, whoever
277! DSAVG assumes IFFLOW to work correctly. IFHEAT still doesn't work with CMT-nek
278 include 'SIZE'
279 include 'TOTAL'
280 include 'CMTDATA'
281 include 'DG'
282 integer lfq,heresize,hdsize
283 parameter (lxyz=lx1*ly1*lz1)
284 real residual(lxyz,nelt)
285 integer e,f
286
287 nxyz =lx1*ly1*lz1
288 nxz =lx1*lz1
289 nface=2*ldim
290 nxzf =nxz*nface
291 nfq =nxzf*nelt
292
293! ensure continuity at faces. doing this before |abs| causes some cancellation that,
294! so far, appears to beneficially reduce spikiness at faces.
295
296 do e=1,nelt
297 do i=1,nxyz
298 residual(i,e)=residual(i,e)*meshh(e)**2
299 enddo
300 enddo
301
302 call dsavg(residual) ! signed, can cancel at faces. Hope it does
303
304 do e=1,nelt
305 do i=1,nxyz
306 residual(i,e)=abs(residual(i,e))
307 enddo
308 enddo
309
310 call cmult(residual,c_sub_e,nxyz*nelt)
311
312 if (maxdiff .ne. 0) then
313 const=1.0/maxdiff
314 call cmult(residual,const,nxyz*nelt)
315 endif
316
317 return
318 end
319
320!-----------------------------------------------------------------------
321
322 subroutine evmsmooth(resvisc,wavevisc,endpoints)
323 include 'SIZE'
324 include 'INPUT'
325 real resvisc(lx1,ly1,lz1,nelt),wavevisc(lx1,ly1,lz1,nelt)
326 real rtmp
327 common /ctmp1/ rtmp(lx1,ly1,lz1)
328! are faces included in smoothing? if not (say they're fixed by dsavg) then
329! endpoints should be .false.
330 logical endpoints
331 integer e
332! clip residual viscosity and smooth it.
333! actually just smooth first and then clip it. Just average-with-my-neighbors for now
334! weight the neighbors according to GLL-spacing instead of uniform-grid
335! at a later date.
336
337 nxyz=lx1*ly1*lz1
338 kstart=1
339 kend=lz1
340 rldim=1.0/ldim
341
342 if (endpoints) then
343 istart=1
344 jstart=1
345 iend=lx1
346 jend=ly1
347 else
348 istart=2
349 jstart=2
350 iend=lx1-1
351 jend=ly1-1
352 if (if3d) then
353 kstart=2
354 kend=lz1-1
355 endif
356 endif
357
358! yes I know this loop sucks and I need to write a matrix or something
359 do e=1,nelt
360 do i=1,nxyz
361 if (wavevisc(i,1,1,e).le. resvisc(i,1,1,e))
362 > resvisc(i,1,1,e)= wavevisc(i,1,1,e)
363 enddo
364 call copy(rtmp,resvisc(1,1,1,e),nxyz) ! really only for .false.
365 do iz=kstart,kend
366 if (if3d) then
367 km1=iz-1
368 kp1=iz+1
369 izm=km1
370! if (km1 .lt. 1) izm=iz ! bias towards {{face point}}
371 if (km1 .lt. 1) izm=kp1 ! Guermond symmetry
372 izp=kp1
373! if (kp1 .gt. lz1) izp=iz ! bias towards {{face point}}
374 if (kp1 .gt. lz1) izp=km1 ! Guermond symmetry
375 else
376 izm=iz
377 izp=iz
378 endif
379 do iy=jstart,jend
380 jm1=iy-1
381 jp1=iy+1
382 iym=jm1
383! if (jm1 .lt. 1) iym=iy ! bias towards {{face point}}
384 if (jm1 .lt. 1) iym=jp1 ! Guermond symmetry
385 iyp=jp1
386! if (jp1 .gt. ly1) iyp=iy ! bias toward {{face point}}
387 if (jp1 .gt. ly1) iyp=jm1 ! Guermond symmetry
388 do ix=istart,iend
389 im1=ix-1
390 ip1=ix+1
391 ixm=im1
392! if (im1 .lt. 1) ixm=ix ! bias towards {{face point}}
393 if (im1 .lt. 1) ixm=ip1 ! Guermond symmetry
394 ixp=ip1
395! if (ip1 .gt. lx1) ixp=ix ! bias towards {{face point}}
396 if (ip1 .gt. lx1) ixp=im1 ! Guermond symmetry
397 x0 = resvisc(ix ,iy ,iz ,e)
398 x1 = resvisc(ixm,iy ,iz ,e)
399 x2 = resvisc(ixp,iy ,iz ,e)
400 x3 = resvisc(ix ,iym,iz ,e)
401 x4 = resvisc(ix ,iyp,iz ,e)
402 if (if3d) then
403 x5 = resvisc(ix ,iy ,izm,e)
404 x6 = resvisc(ix ,iy ,izp,e)
405 else
406 x5=0.0
407 x6=0.0
408 endif
409 rtmp(ix,iy,iz)=0.25*(2.0*ldim*x0+x1+x2+x3+x4+x5+x6)
410 > *rldim
411 enddo
412 enddo
413 enddo
414 call copy(resvisc(1,1,1,e),rtmp,nxyz)
415 enddo
416
417 return
418 end
419
420!-----------------------------------------------------------------------
421
422 subroutine wavevisc(numax)
423! max entropy visc, defined by Guermond, Popov, whoever (chapter-verse)
424! as
425! numax = c_max*h*max(dH/dU)
426! in a given element
427 include 'SIZE'
428 include 'TOTAL'
429 include 'CMTDATA'
430 parameter (lxyz=lx1*ly1*lz1)
431 common /scrns/ wavespeed(lxyz)
432 real wavespeed
433 real maxeig
434 real numax(lxyz,nelt)
435 integer e
436
437 nxyz=lx1*ly1*lz1
438
439 do e=1,nelt
440 do i=1,nxyz
441 wavespeed(i)=csound(i,1,1,e)+
442 > sqrt(vx(i,1,1,e)**2+vy(i,1,1,e)**2+vz(i,1,1,e)**2)
443 enddo
444 maxeig=vlamax(wavespeed,nxyz)
445! Zingan (2015) only. not long for this world
446! rhomax(e)=vlamax(vtrans(1,1,1,e,irho),nxyz)
447 do i=1,nxyz
448 numax(i,e)=c_max*maxeig*meshh(e)
449 enddo
450 enddo
451
452 call max_to_trilin(t(1,1,1,1,3))
453
454 return
455 end
456
457!-----------------------------------------------------------------------
458
459 subroutine max_to_trilin(field)
460! stupid subroutine to take a stupid uniform field and compute a trilinear
461! tent between maximum shared values at the vertices.
462 include 'SIZE'
463 include 'TOTAL'
464 real field(lx1,ly1,lz1,nelt)
465 integer e
466
467 nxyz=lx1*ly1*lz1
468
469! get maxima on faces
470 call dsop(field,'MAX',lx1,ly1,lz1)
471
472! trilinear interpolation. you should adapt xyzlin to your needs instead
473 do e=1,nelt
474 p000=field(1, 1, 1, e)
475 p100=field(lx1,1, 1, e)
476 p010=field(1, ly1,1, e)
477 p110=field(lx1,ly1,1, e)
478 p001=field(1, 1, lz1,e)
479 p101=field(lx1,1, lz1,e)
480 p011=field(1, ly1,lz1,e)
481 p111=field(lx1,ly1,lz1,e)
482 c1=p100-p000
483 c2=p010-p000
484 c3=p001-p000
485 c4=p110-p010-p100+p000
486 c5=p011-p001-p010+p000
487 c6=p101-p001-p100+p000
488 c7=p111-p011-p101-p110+p100+p001+p010-p000
489 rdx=1.0/(xm1(lx1,1,1,e)-xm1(1,1,1,e)) ! cubes only!!!
490 rdy=1.0/(ym1(1,ly1,1,e)-ym1(1,1,1,e))
491 rdz=0.0
492 if(if3d) rdz=1.0/(zm1(1,1,lz1,e)-zm1(1,1,1,e))
493 do i=1,nxyz
494 deltax=rdx*(xm1(i,1,1,e)-xm1(1,1,1,e)) ! cubes only!!!
495 deltay=rdy*(ym1(i,1,1,e)-ym1(1,1,1,e))
496 deltaz=0.0
497 if (if3d) deltaz=rdz*(zm1(i,1,1,e)-zm1(1,1,1,e))
498 field(i,1,1,e)=p000+c1*deltax+c2*deltay+c3*deltaz+
499 > c4*deltax*deltay+c5*deltay*deltaz+
500 > c6*deltaz*deltax+c7*deltay*deltaz*deltax
501 enddo
502 enddo
503
504 return
505 end
506
507!-----------------------------------------------------------------------
508
509 subroutine perssonperaire(shkdet,var,shtmp)
510! Peraire & Persson (2006) Eq. 7
511 include 'SIZE'
512 include 'TOTAL'
513 include 'CMTDATA'
514 parameter (lxyz=lx1*ly1*lz1)
515 common /scrns/ fvar (lxyz,lelt)
516! $ , ytm1 (lx1,ly1,lz1,lelv)
517 real shkdet(lxyz,nelt),var(lxyz,*),shtmp(lxyz,*)
518 integer e
519! this whole common block is atrocious
520 common /CMTFILTERS/ intv(lx1,lx1),intt(lx1,lx1)
521 $ , intvd(lxd,lxd),inttd(lxd,lxd)
522 $ , wk1(lx1,lx1,lz1),wk2(lx1,lx1,lz1)
523 $ , wkd1(lxd,lxd,lzd),wkd2(lxd,lxd,lzd)
524 real intv, intt, intvd, inttd
525 real kappa
526 integer icalld
527 save icalld
528 data icalld/0/
529
530 if (icalld.eq.0) then
531 ncut=-110
532 icalld=1
533 call userfilt(intv,zgm1,lx1,nio)
534 endif
535
536 n=lx1*ly1*lz1
537 ntot=n*nelt
538
539! JH082418 Hardcoded threshold
540 cpp=0.001
541
542 call copy(fvar,var,ntot)
543 call filterq(fvar,intv,lx1,lz1,wk1,wk2,intt,if3d,dmax)
544
545 do e=1,nelt
546 call sub3(shkdet(1,e),var(1,e),fvar(1,e),n) ! store u-u^ in shkdet
547 shock_or_not=vlsc3(shkdet(1,e),shkdet(1,e),bm1(1,1,1,e),n)
548 denom=vlsc3(var(1,e),var(1,e),bm1(1,1,1,e),n)
549 if (denom .gt. 0.0) then ! mask this somehow. ifs are bad
550 shock_or_not=(shock_or_not/denom)*(lx1-1)**4
551! shock_or_not=log10(shock_or_not)
552 if (shock_or_not .gt. cpp) then
553! call cfill(shkdet(1,e),shock_or_not,n)
554 call rone(shkdet(1,e),n)
555 else
556 call rzero(shkdet(1,e),n)
557 endif
558 else
559 write(6,*) 'nid,iel=',nid,e
560 write(6,*) 'variables that change sign are bad for shox'
561 call exitt
562 endif
563 enddo
564
565 return
566 end
Note: See TracBrowser for help on using the repository browser.