1 |
C $Header: $ |
2 |
C $Name: $ |
3 |
|
4 |
# include "GAD_OPTIONS.h" |
5 |
|
6 |
SUBROUTINE GAD_OSC_MUL_R(ir,hh,mask,ohat,scal) |
7 |
C |================================================================| |
8 |
C | OSC_MUL_R: evaluate WENO oscillation weights in R. | |
9 |
C |================================================================| |
10 |
|
11 |
implicit none |
12 |
|
13 |
C =============================================== global variables |
14 |
# include "SIZE.h" |
15 |
|
16 |
integer ir,hh |
17 |
_RL mask(1-3:Nr+3) |
18 |
_RL ohat(1:2, |
19 |
& 1-3:Nr+3) |
20 |
_RL scal(1:2) |
21 |
|
22 |
integer ii |
23 |
_RL dels,dfs1,dfs2 |
24 |
_RL osum,zero,mval |
25 |
_RL oval,omin,omax |
26 |
|
27 |
C =============================== calc. WENO oscillation weighting |
28 |
c omin = +huge(+1. _d 0) |
29 |
c omax = -huge(+1. _d 0) |
30 |
omin = +1. _d 99 |
31 |
omax = -1. _d 99 |
32 |
|
33 |
zero = 1. _d -20 |
34 |
mval = 1. _d + 0 |
35 |
|
36 |
do ii = ir-hh, ir+hh |
37 |
|
38 |
C =============================== calc. derivatives centred on II. |
39 |
dels = (ii - ir) * 2. _d 0 |
40 |
|
41 |
dfs1 = ohat(1,ii) |
42 |
dfs2 = ohat(2,ii) |
43 |
|
44 |
dfs1 = dfs1 + dfs2 * dels |
45 |
|
46 |
C =============================== oscl. = NORM(H^N * D^N/DR^N(F)). |
47 |
oval = (2. _d 0 * dfs1)**2 |
48 |
& + (4. _d 0 * dfs2)**2 |
49 |
|
50 |
if (oval.lt.omin) omin = oval |
51 |
if (oval.gt.omax) omax = oval |
52 |
|
53 |
C =============================== any mask across oscil. stencil |
54 |
mval = mval * mask(ii) |
55 |
|
56 |
end do |
57 |
|
58 |
if (mval .gt. 0. _d 0) then |
59 |
|
60 |
C =============================== calc. WENO-style profile weights |
61 |
scal(1) = 1. _d 5 |
62 |
& / (omax + zero)**3 |
63 |
scal(2) = 1. _d 0 |
64 |
& / (omin + zero)**3 |
65 |
|
66 |
osum = scal(1) + scal(2) |
67 |
scal(1) = scal(1) / osum |
68 |
scal(2) = scal(2) / osum |
69 |
|
70 |
else |
71 |
|
72 |
C =============================== default to MONO. profile weights |
73 |
scal(1) = 0. _d 0 |
74 |
scal(2) = 1. _d 0 |
75 |
|
76 |
end if |
77 |
|
78 |
return |
79 |
|
80 |
c end subroutine GAD_OSC_MUL_R |
81 |
end |