1 |
C $Header: $ |
2 |
C $Name: $ |
3 |
|
4 |
# include "GAD_OPTIONS.h" |
5 |
|
6 |
SUBROUTINE GAD_PQM_P5E_Y(bi,bj,kk,ix, |
7 |
& mask,fbar,edge,myThid) |
8 |
C |================================================================| |
9 |
C | PQM_P5E_Y: approximate edge values with degree-5 polynomials. | |
10 |
C | Fixed grid-spacing variant in Y. | |
11 |
C |================================================================| |
12 |
|
13 |
implicit none |
14 |
|
15 |
C =============================================== global variables |
16 |
# include "SIZE.h" |
17 |
# include "GRID.h" |
18 |
# include "GAD.h" |
19 |
|
20 |
C ====================================================== arguments |
21 |
integer bi,bj,kk,ix |
22 |
_RL mask(1-OLy:sNy+OLy) |
23 |
_RL fbar(1-OLy:sNy+OLy) |
24 |
_RL edge(1:2, |
25 |
& 1-OLy:sNy+OLy) |
26 |
integer myThid |
27 |
|
28 |
C ====================================================== variables |
29 |
integer iy |
30 |
_RL mloc(-3:+2) |
31 |
_RL floc(-3:+2) |
32 |
_RL ftmp |
33 |
|
34 |
C ==================== reconstruct 5th--order accurate edge values |
35 |
do iy = 1-OLy+3, sNy+OLy-2 |
36 |
|
37 |
C ================ mask local stencil: expand from centre outwards |
38 |
mloc(-1) = mask(iy-1) |
39 |
mloc(+0) = mask(iy+0) |
40 |
|
41 |
floc(-1) = fbar(iy+0) |
42 |
& + mloc(-1)*(fbar(iy-1)-fbar(iy+0)) |
43 |
floc(+0) = fbar(iy-1) |
44 |
& + mloc(+0)*(fbar(iy+0)-fbar(iy-1)) |
45 |
|
46 |
mloc(-2) = mask(iy-2) * mloc(-1) |
47 |
mloc(-3) = mask(iy-3) * mloc(-2) |
48 |
|
49 |
ftmp = 2. _d 0 * floc(-1) - floc(+0) |
50 |
floc(-2) = ftmp |
51 |
& + mloc(-2)*(fbar(iy-2)-ftmp) |
52 |
ftmp = 2. _d 0 * floc(-2) - floc(-1) |
53 |
floc(-3) = ftmp |
54 |
& + mloc(-3)*(fbar(iy-3)-ftmp) |
55 |
|
56 |
mloc(+1) = mask(iy+1) * mloc(+0) |
57 |
mloc(+2) = mask(iy+2) * mloc(+1) |
58 |
|
59 |
ftmp = 2. _d 0 * floc(+0) - floc(-1) |
60 |
floc(+1) = ftmp |
61 |
& + mloc(+1)*(fbar(iy+1)-ftmp) |
62 |
ftmp = 2. _d 0 * floc(+1) - floc(+0) |
63 |
floc(+2) = ftmp |
64 |
& + mloc(+2)*(fbar(iy+2)-ftmp) |
65 |
|
66 |
C ================ centred, 5th-order interpolation for edge value |
67 |
edge(1,iy) = |
68 |
& +( 1. _d 0/60. _d 0)*(floc(-3)+floc(+2)) |
69 |
& -( 8. _d 0/60. _d 0)*(floc(-2)+floc(+1)) |
70 |
& +(37. _d 0/60. _d 0)*(floc(-1)+floc(+0)) |
71 |
|
72 |
edge(2,iy) = ( |
73 |
& -( 1. _d 0/90. _d 0)*(floc(-3)-floc(+2)) |
74 |
& +( 5. _d 0/36. _d 0)*(floc(-2)-floc(+1)) |
75 |
& -(49. _d 0/36. _d 0)*(floc(-1)-floc(+0)) |
76 |
& ) * recip_dyC(ix,iy,bi,bj) |
77 |
|
78 |
end do |
79 |
|
80 |
return |
81 |
|
82 |
c end subroutine GAD_PQM_P5E_Y |
83 |
end |