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

Annotation of /MITgcm/pkg/generic_advdiff/gad_ppm_p3e_x.F

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


Revision 1.1 - (hide annotations) (download)
Sun Mar 13 01:44:03 2016 UTC (8 years, 2 months 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 jmc 1.1 C $Header: $
2     C $Name: $
3    
4     # include "GAD_OPTIONS.h"
5    
6     SUBROUTINE GAD_PPM_P3E_X(bi,bj,kk,iy,
7     & mask,fbar,edge,myThid)
8     C |================================================================|
9     C | PPM_P3E_X: approximate edge values with degree-3 polynomials. |
10     C | Fixed grid-spacing variant in X. |
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,iy
22     _RL mask(1-OLx:sNx+OLx)
23     _RL fbar(1-OLx:sNx+OLx)
24     _RL edge(1-OLx:sNx+OLx)
25     integer myThid
26    
27     C ====================================================== variables
28     integer ix
29     _RL mloc(-2:+1)
30     _RL floc(-2:+1)
31     _RL ftmp
32    
33     do ix = 1-OLx+2, sNx+OLx-1
34    
35     C ================ mask local stencil: expand from centre outwards
36     mloc(-1) = mask(ix-1)
37     mloc(+0) = mask(ix+0)
38    
39     floc(-1) = fbar(ix+0)
40     & + mloc(-1)*(fbar(ix-1)-fbar(ix+0))
41     floc(+0) = fbar(ix-1)
42     & + mloc(+0)*(fbar(ix+0)-fbar(ix-1))
43    
44     mloc(-2) = mask(ix-2) * mloc(-1)
45    
46     ftmp = 2. _d 0 * floc(-1) - floc(+0)
47     floc(-2) = ftmp
48     & + mloc(-2)*(fbar(ix-2)-ftmp)
49    
50     mloc(+1) = mask(ix+1) * mloc(+0)
51    
52     ftmp = 2. _d 0 * floc(+0) - floc(-1)
53     floc(+1) = ftmp
54     & + mloc(+1)*(fbar(ix+1)-ftmp)
55    
56     C ================ centred, 3rd-order interpolation for edge value
57     edge(ix) =
58     & -(1. _d 0 / 12. _d 0)*(floc(-2)+floc(+1))
59     & +(7. _d 0 / 12. _d 0)*(floc(-1)+floc(+0))
60    
61     end do
62    
63     return
64    
65     c end subroutine GAD_PPM_P3E_X
66     end

  ViewVC Help
Powered by ViewVC 1.1.22