/[MITgcm]/MITgcm/pkg/generic_advdiff/gad_ppm_p3e_y.F
ViewVC logotype

Contents of /MITgcm/pkg/generic_advdiff/gad_ppm_p3e_y.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.1 - (show annotations) (download)
Sun Mar 13 01:44:03 2016 UTC (8 years, 1 month ago) by jmc
Branch: MAIN
CVS Tags: checkpoint66g, checkpoint66f, checkpoint66e, checkpoint66d, checkpoint66c, checkpoint66b, checkpoint66a, checkpoint66o, checkpoint66n, checkpoint66m, checkpoint66l, checkpoint66k, checkpoint66j, checkpoint66i, checkpoint66h, checkpoint65z, checkpoint65x, checkpoint65y, checkpoint65v, checkpoint65w, checkpoint65u, HEAD
- from Darren: add PPM and PQM advection schemes (number 40-42 and 50-52)
  with 2 types of limiter (see: Engwirda & Kelley, submit. to JCP);
  Note (from Darren): unlimited PPM/PQM scheme (40 & 50) are just for
  testing and not for actual use.

1 C $Header: $
2 C $Name: $
3
4 # include "GAD_OPTIONS.h"
5
6 SUBROUTINE GAD_PPM_P3E_Y(bi,bj,kk,ix,
7 & mask,fbar,edge,myThid)
8 C |================================================================|
9 C | PPM_P3E_Y: approximate edge values with degree-3 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-OLy:sNy+OLy)
25 integer myThid
26
27 C ====================================================== variables
28 integer iy
29 _RL mloc(-2:+1)
30 _RL floc(-2:+1)
31 _RL ftmp
32
33 C ==================== reconstruct 3rd--order accurate edge values
34 do iy = 1-OLy+2, sNy+OLy-1
35
36 C ================ mask local stencil: expand from centre outwards
37 mloc(-1) = mask(iy-1)
38 mloc(+0) = mask(iy+0)
39
40 floc(-1) = fbar(iy+0)
41 & + mloc(-1)*(fbar(iy-1)-fbar(iy+0))
42 floc(+0) = fbar(iy-1)
43 & + mloc(+0)*(fbar(iy+0)-fbar(iy-1))
44
45 mloc(-2) = mask(iy-2) * mloc(-1)
46
47 ftmp = 2. _d 0 * floc(-1) - floc(+0)
48 floc(-2) = ftmp
49 & + mloc(-2)*(fbar(iy-2)-ftmp)
50
51 mloc(+1) = mask(iy+1) * mloc(+0)
52
53 ftmp = 2. _d 0 * floc(+0) - floc(-1)
54 floc(+1) = ftmp
55 & + mloc(+1)*(fbar(iy+1)-ftmp)
56
57 C ================ centred, 5th-order interpolation for edge value
58 edge(iy) =
59 & -(1. _d 0 / 12. _d 0)*(floc(-2)+floc(+1))
60 & +(7. _d 0 / 12. _d 0)*(floc(-1)+floc(+0))
61
62 end do
63
64 return
65
66 c end subroutine GAD_PPM_P3E_Y
67 end

  ViewVC Help
Powered by ViewVC 1.1.22