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

Annotation of /MITgcm/pkg/generic_advdiff/gad_pqm_p5e_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, 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 jmc 1.1 C $Header: $
2     C $Name: $
3    
4     # include "GAD_OPTIONS.h"
5    
6     SUBROUTINE GAD_PQM_P5E_X(bi,bj,kk,iy,
7     & mask,fbar,edge,myThid)
8     C |================================================================|
9     C | PQM_P5E_X: approximate edge values with degree-5 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:2,
25     & 1-OLx:sNx+OLx)
26     integer myThid
27    
28     C ====================================================== variables
29     integer ix
30     _RL mloc(-3:+2)
31     _RL floc(-3:+2)
32     _RL ftmp
33    
34     do ix = 1-OLx+3, sNx+OLx-2
35    
36     C ================ mask local stencil: expand from centre outwards
37     mloc(-1) = mask(ix-1)
38     mloc(+0) = mask(ix+0)
39    
40     floc(-1) = fbar(ix+0)
41     & + mloc(-1)*(fbar(ix-1)-fbar(ix+0))
42     floc(+0) = fbar(ix-1)
43     & + mloc(+0)*(fbar(ix+0)-fbar(ix-1))
44    
45     mloc(-2) = mask(ix-2) * mloc(-1)
46     mloc(-3) = mask(ix-3) * mloc(-2)
47    
48     ftmp = 2. _d 0 * floc(-1) - floc(+0)
49     floc(-2) = ftmp
50     & + mloc(-2)*(fbar(ix-2)-ftmp)
51     ftmp = 2. _d 0 * floc(-2) - floc(-1)
52     floc(-3) = ftmp
53     & + mloc(-3)*(fbar(ix-3)-ftmp)
54    
55     mloc(+1) = mask(ix+1) * mloc(+0)
56     mloc(+2) = mask(ix+2) * mloc(+1)
57    
58     ftmp = 2. _d 0 * floc(+0) - floc(-1)
59     floc(+1) = ftmp
60     & + mloc(+1)*(fbar(ix+1)-ftmp)
61     ftmp = 2. _d 0 * floc(+1) - floc(+0)
62     floc(+2) = ftmp
63     & + mloc(+2)*(fbar(ix+2)-ftmp)
64    
65     C ================ centred, 5th-order interpolation for edge value
66     edge(1,ix) =
67     & +( 1. _d 0/60. _d 0)*(floc(-3)+floc(+2))
68     & -( 8. _d 0/60. _d 0)*(floc(-2)+floc(+1))
69     & +(37. _d 0/60. _d 0)*(floc(-1)+floc(+0))
70    
71     edge(2,ix) = (
72     & -( 1. _d 0/90. _d 0)*(floc(-3)-floc(+2))
73     & +( 5. _d 0/36. _d 0)*(floc(-2)-floc(+1))
74     & -(49. _d 0/36. _d 0)*(floc(-1)-floc(+0))
75     & ) * recip_dxC(ix,iy,bi,bj)
76    
77     end do
78    
79     return
80    
81     c end subroutine GAD_PQM_P5E_X
82     end

  ViewVC Help
Powered by ViewVC 1.1.22