MOM6
plm_functions Module Reference

Functions/Subroutines

subroutine, public plm_reconstruction (N, h, u, ppoly_E, ppoly_coefficients)
 
subroutine, public plm_boundary_extrapolation (N, h, u, ppoly_E, ppoly_coefficients)
 

Variables

real, parameter h_neglect = 1.E-30
 

Function/Subroutine Documentation

◆ plm_boundary_extrapolation()

subroutine, public plm_functions::plm_boundary_extrapolation ( integer, intent(in)  N,
real, dimension(:), intent(in)  h,
real, dimension(:), intent(in)  u,
real, dimension(:,:), intent(inout)  ppoly_E,
real, dimension(:,:), intent(inout)  ppoly_coefficients 
)

Definition at line 212 of file PLM_functions.F90.

References h_neglect.

Referenced by mom_remapping::build_reconstructions_1d(), mom_ale::pressure_gradient_plm(), and regrid_interp::regridding_set_ppolys().

212 !------------------------------------------------------------------------------
213 ! Reconstruction by linear polynomials within boundary cells.
214 ! The left and right edge values in the left and right boundary cells,
215 ! respectively, are estimated using a linear extrapolation within the cells.
216 !
217 ! This extrapolation is EXACT when the underlying profile is linear.
218 !
219 ! N: number of cells in grid
220 ! h: thicknesses of grid cells
221 ! u: cell averages to use in constructing piecewise polynomials
222 ! ppoly_E : edge values of piecewise polynomials
223 ! ppoly_coefficients : coefficients of piecewise polynomials
224 !
225 ! It is assumed that the size of the array 'u' is equal to the number of cells
226 ! defining 'grid' and 'ppoly'. No consistency check is performed here.
227 !------------------------------------------------------------------------------
228 
229  ! Arguments
230  integer, intent(in) :: n ! Number of cells
231  real, dimension(:), intent(in) :: h ! cell widths (size N)
232  real, dimension(:), intent(in) :: u ! cell averages (size N)
233  real, dimension(:,:), intent(inout) :: ppoly_e
234  real, dimension(:,:), intent(inout) :: ppoly_coefficients
235 
236  ! Local variables
237  real :: u0, u1 ! cell averages
238  real :: h0, h1 ! corresponding cell widths
239  real :: slope ! retained PLM slope
240 
241  ! -----------------------------------------
242  ! Left edge value in the left boundary cell
243  ! -----------------------------------------
244  h0 = h(1) + h_neglect
245  h1 = h(2) + h_neglect
246 
247  u0 = u(1)
248  u1 = u(2)
249 
250  ! The h2 scheme is used to compute the right edge value
251  ppoly_e(1,2) = (u0*h1 + u1*h0) / (h0 + h1)
252 
253  ! The standard PLM slope is computed as a first estimate for the
254  ! reconstruction within the cell
255  slope = 2.0 * ( ppoly_e(1,2) - u0 )
256 
257  ppoly_e(1,1) = u0 - 0.5 * slope
258  ppoly_e(1,2) = u0 + 0.5 * slope
259 
260  ppoly_coefficients(1,1) = ppoly_e(1,1)
261  ppoly_coefficients(1,2) = ppoly_e(1,2) - ppoly_e(1,1)
262 
263  ! ------------------------------------------
264  ! Right edge value in the left boundary cell
265  ! ------------------------------------------
266  h0 = h(n-1) + h_neglect
267  h1 = h(n) + h_neglect
268 
269  u0 = u(n-1)
270  u1 = u(n)
271 
272  ! The h2 scheme is used to compute the right edge value
273  ppoly_e(n,1) = (u0*h1 + u1*h0) / (h0 + h1)
274 
275  ! The standard PLM slope is computed as a first estimate for the
276  ! reconstruction within the cell
277  slope = 2.0 * ( u1 - ppoly_e(n,1) )
278 
279  ppoly_e(n,1) = u1 - 0.5 * slope
280  ppoly_e(n,2) = u1 + 0.5 * slope
281 
282  ppoly_coefficients(n,1) = ppoly_e(n,1)
283  ppoly_coefficients(n,2) = ppoly_e(n,2) - ppoly_e(n,1)
284 
Here is the caller graph for this function:

◆ plm_reconstruction()

subroutine, public plm_functions::plm_reconstruction ( integer, intent(in)  N,
real, dimension(:), intent(in)  h,
real, dimension(:), intent(in)  u,
real, dimension(:,:), intent(inout)  ppoly_E,
real, dimension(:,:), intent(inout)  ppoly_coefficients 
)

Definition at line 26 of file PLM_functions.F90.

References h_neglect.

Referenced by mom_remapping::build_reconstructions_1d(), mom_ale::pressure_gradient_plm(), regrid_interp::regridding_set_ppolys(), and mom_remapping::remapping_unit_tests().

26 !------------------------------------------------------------------------------
27 ! Reconstruction by linear polynomials within each cell.
28 !
29 ! N: number of cells in grid
30 ! h: thicknesses of grid cells
31 ! u: cell averages to use in constructing piecewise polynomials
32 ! ppoly_E : edge values of piecewise polynomials
33 ! ppoly_coefficients : coefficients of piecewise polynomials
34 !
35 ! It is assumed that the size of the array 'u' is equal to the number of cells
36 ! defining 'grid' and 'ppoly'. No consistency check is performed here.
37 !------------------------------------------------------------------------------
38 
39  ! Arguments
40  integer, intent(in) :: n ! Number of cells
41  real, dimension(:), intent(in) :: h ! cell widths (size N)
42  real, dimension(:), intent(in) :: u ! cell averages (size N)
43  real, dimension(:,:), intent(inout) :: ppoly_e
44  real, dimension(:,:), intent(inout) :: ppoly_coefficients
45 
46  ! Local variables
47  integer :: k ! loop index
48  real :: u_l, u_c, u_r ! left, center and right cell averages
49  real :: h_l, h_c, h_r, h_cn ! left, center and right cell widths
50  real :: sigma_l, sigma_c, sigma_r ! left, center and right
51  ! van Leer slopes
52  real :: slope ! retained PLM slope
53  real :: a, b ! auxiliary variables
54  real :: u_min, u_max, e_l, e_r, edge
55  real :: almost_one, almost_two
56  real, dimension(N) :: slp, mslp
57 
58  almost_one = 1. - epsilon(slope)
59  almost_two = 2. * almost_one
60 
61  ! Loop on interior cells
62  do k = 2,n-1
63 
64  ! Get cell averages
65  u_l = u(k-1) ; u_c = u(k) ; u_r = u(k+1)
66 
67  ! Get cell widths
68  h_l = h(k-1) ; h_c = h(k) ; h_r = h(k+1)
69  h_cn = max( h_c, h_neglect ) ! To avoid division by zero
70 
71  ! Side differences
72  sigma_r = u_r - u_c
73  sigma_l = u_c - u_l
74 
75  ! This is the second order slope given by equation 1.7 of
76  ! Piecewise Parabolic Method, Colella and Woodward (1984),
77  ! http://dx.doi.org/10.1016/0021-991(84)90143-8.
78  ! For uniform resolution it simplifies to ( u_r - u_l )/2 .
79  ! sigma_c = ( h_c / ( h_cn + ( h_l + h_r ) ) ) * ( &
80  ! ( 2.*h_l + h_c ) / ( h_r + h_cn ) * sigma_r &
81  ! + ( 2.*h_r + h_c ) / ( h_l + h_cn ) * sigma_l )
82 
83  ! This is the original estimate of the second order slope from Laurent
84  ! but multiplied by h_c
85  sigma_c = 2.0 * ( u_r - u_l ) * ( h_c / ( h_l + 2.0*h_c + h_r + h_neglect) )
86 
87  if ( (sigma_l * sigma_r) > 0.0 ) then
88  ! This limits the slope so that the edge values are bounded by the
89  ! two cell averages spanning the edge.
90  u_min = min( u_l, u_c, u_r )
91  u_max = max( u_l, u_c, u_r )
92  slope = sign( min( abs(sigma_c), 2.*min( u_c - u_min, u_max - u_c ) ), sigma_c )
93  else
94  ! Extrema in the mean values require a PCM reconstruction avoid generating
95  ! larger extreme values.
96  slope = 0.0
97  end if
98 
99  ! This block tests to see if roundoff causes edge values to be out of bounds
100  u_min = min( u_l, u_c, u_r )
101  u_max = max( u_l, u_c, u_r )
102  if (u_c - 0.5*abs(slope) < u_min .or. u_c + 0.5*abs(slope) > u_max) then
103  slope = slope * almost_one
104  endif
105 
106  ! An attempt to avoid inconsistency when the values become unrepresentable.
107  if (abs(slope) < 1.e-140) slope = 0.
108 
109  ! Safety check - this block really should not be needed ...
110 ! if (u_c - 0.5*abs(slope) < u_min .or. u_c + 0.5*abs(slope) > u_max) then
111 ! write(0,*) 'l,c,r=',u_l,u_c,u_r
112 ! write(0,*) 'min,max=',u_min,u_max
113 ! write(0,*) 'slp=',slope
114 ! sigma_l = u_c-0.5*abs(slope)
115 ! sigma_r = u_c+0.5*abs(slope)
116 ! write(0,*) 'lo,hi=',sigma_l,sigma_r
117 ! write(0,*) 'elo,ehi=',sigma_l-u_min,sigma_r-u_max
118 ! stop 'Limiter failed!'
119 ! endif
120 
121  slp(k) = slope
122  ppoly_e(k,1) = u_c - 0.5 * slope
123  ppoly_e(k,2) = u_c + 0.5 * slope
124 
125  end do ! end loop on interior cells
126 
127  ! Boundary cells use PCM. Extrapolation is handled in a later routine.
128  slp(1) = 0.
129  ppoly_e(1,2) = u(1)
130  slp(n) = 0.
131  ppoly_e(n,1) = u(n)
132 
133  ! This loop adjusts the slope so that edge values are monotonic.
134  do k = 2, n-1
135  u_l = u(k-1) ; u_c = u(k) ; u_r = u(k+1)
136  e_r = ppoly_e(k-1,2) ! Right edge from cell k-1
137  e_l = ppoly_e(k+1,1) ! Left edge from cell k
138  mslp(k) = abs(slp(k))
139  u_min = min(e_r, u_c)
140  u_max = max(e_r, u_c)
141  edge = u_c - 0.5 * slp(k)
142  if ( ( edge - e_r ) * ( u_c - edge ) < 0. ) then
143  edge = 0.5 * ( edge + e_r ) ! * almost_one?
144  mslp(k) = min( mslp(k), abs( edge - u_c ) * almost_two )
145  endif
146  edge = u_c + 0.5 * slp(k)
147  if ( ( edge - u_c ) * ( e_l - edge ) < 0. ) then
148  edge = 0.5 * ( edge + e_l ) ! * almost_one?
149  mslp(k) = min( mslp(k), abs( edge - u_c ) * almost_two )
150  endif
151  enddo ! end loop on interior cells
152  mslp(1) = 0.
153  mslp(n) = 0.
154 
155  ! Check that the above adjustment worked
156 ! do K = 2, N-1
157 ! u_r = u(k-1) + 0.5 * sign( mslp(k-1), slp(k-1) ) ! Right edge from cell k-1
158 ! u_l = u(k) - 0.5 * sign( mslp(k), slp(k) ) ! Left edge from cell k
159 ! if ( (u(k)-u(k-1)) * (u_l-u_r) < 0. ) then
160 ! stop 'Adjustment failed!'
161 ! endif
162 ! enddo ! end loop on interior cells
163 
164  ! Store and return edge values and polynomial coefficients.
165  ppoly_e(1,1) = u(1)
166  ppoly_e(1,2) = u(1)
167  ppoly_coefficients(1,1) = u(1)
168  ppoly_coefficients(1,2) = 0.
169  do k = 2, n-1
170  slope = sign( mslp(k), slp(k) )
171  u_l = u(k) - 0.5 * slope ! Left edge value of cell k
172  u_r = u(k) + 0.5 * slope ! Right edge value of cell k
173 
174  ! Check that final edge values are bounded
175  u_min = min( u(k-1), u(k) )
176  u_max = max( u(k-1), u(k) )
177  if (u_l<u_min .or. u_l>u_max) then
178  write(0,*) 'u(k-1)=',u(k-1),'u(k)=',u(k),'slp=',slp(k),'u_l=',u_l
179  stop 'Left edge out of bounds'
180  endif
181  u_min = min( u(k+1), u(k) )
182  u_max = max( u(k+1), u(k) )
183  if (u_r<u_min .or. u_r>u_max) then
184  write(0,*) 'u(k)=',u(k),'u(k+1)=',u(k+1),'slp=',slp(k),'u_r=',u_r
185  stop 'Right edge out of bounds'
186  endif
187 
188  ppoly_e(k,1) = u_l
189  ppoly_e(k,2) = u_r
190  ppoly_coefficients(k,1) = u_l
191  ppoly_coefficients(k,2) = ( u_r - u_l )
192  ! Check to see if this evaluation of the polynomial at x=1 would be
193  ! monotonic w.r.t. the next cell's edge value. If not, scale back!
194  edge = ppoly_coefficients(k,2) + ppoly_coefficients(k,1)
195  e_r = u(k+1) - 0.5 * sign( mslp(k+1), slp(k+1) )
196  if ( (edge-u(k))*(e_r-edge)<0.) then
197  ppoly_coefficients(k,2) = ppoly_coefficients(k,2) * almost_one
198  endif
199  enddo
200  ppoly_e(n,1) = u(n)
201  ppoly_e(n,2) = u(n)
202  ppoly_coefficients(n,1) = u(n)
203  ppoly_coefficients(n,2) = 0.
204 
Here is the caller graph for this function:

Variable Documentation

◆ h_neglect

real, parameter plm_functions::h_neglect = 1.E-30
private

Definition at line 18 of file PLM_functions.F90.

Referenced by plm_boundary_extrapolation(), and plm_reconstruction().

18 real, parameter :: h_neglect = 1.e-30