Class MUSCLReconstructor
java.lang.Object
neqsim.process.equipment.pipeline.twophasepipe.numerics.MUSCLReconstructor
- All Implemented Interfaces:
Serializable
MUSCL (Monotone Upstream-centered Schemes for Conservation Laws) reconstruction.
Provides second-order spatial accuracy by reconstructing cell-interface values from cell-average values using slope limiters to prevent spurious oscillations near discontinuities.
Reconstruction Formula
For a cell i with value U_i:
- Left interface value: U_{i+1/2,L} = U_i + 0.5 * φ(r) * (U_i - U_{i-1})
- Right interface value: U_{i-1/2,R} = U_i - 0.5 * φ(r) * (U_{i+1} - U_i)
where φ(r) is the slope limiter and r is the ratio of consecutive gradients.
Available Limiters
- Minmod: Most diffusive, very robust
- Van Leer: Good balance of accuracy and stability
- Van Albada: Smooth limiter, good for smooth solutions
- Superbee: Least diffusive, can be oscillatory
- MC (Monotonized Central): Between minmod and superbee
References
- van Leer, B. (1979) - Towards the ultimate conservative difference scheme V.
- Sweby, P.K. (1984) - High resolution schemes using flux limiters
- Version:
- 1.0
- Author:
- Even Solbraa
- See Also:
-
Nested Class Summary
Nested ClassesModifier and TypeClassDescriptionstatic classReconstructed values at a cell interface.static enumAvailable slope limiter types. -
Field Summary
FieldsModifier and TypeFieldDescriptionprivate static final doubleSmall value to prevent division by zero.private MUSCLReconstructor.SlopeLimiterCurrent limiter type.private static final long -
Constructor Summary
ConstructorsConstructorDescriptionDefault constructor using Van Leer limiter.Constructor with specified limiter. -
Method Summary
Modifier and TypeMethodDescriptiondoublecalcLimitedSlope(double deltaLeft, double deltaRight) Calculate limited slope using current limiter.doublecalcLimiter(double r) Calculate slope limiter value φ(r).Get current limiter type.booleanCheck if reconstruction is second-order.doublemc(double r) MC (Monotonized Central) limiter: φ(r) = max(0, min(2, 2r, (1+r)/2)).doubleminmod(double r) Minmod limiter: φ(r) = max(0, min(1, r)).doubleminmod3(double a, double b, double c) Classic minmod function for three arguments.reconstruct(double uLeft, double uCenter, double uRight, double uFarRight) Reconstruct interface values between cells i and i+1.reconstructArray(double[] u) Reconstruct interface values for an array of cell values.voidsetLimiterType(MUSCLReconstructor.SlopeLimiter limiterType) Set limiter type.doublesuperbee(double r) Superbee limiter: φ(r) = max(0, min(2r, 1), min(r, 2)).doublevanAlbada(double r) Van Albada limiter: φ(r) = (r² + r) / (r² + 1).doublevanLeer(double r) Van Leer limiter: φ(r) = (r + |r|) / (1 + |r|).
-
Field Details
-
serialVersionUID
private static final long serialVersionUID- See Also:
-
limiterType
Current limiter type. -
EPSILON
private static final double EPSILONSmall value to prevent division by zero.- See Also:
-
-
Constructor Details
-
MUSCLReconstructor
public MUSCLReconstructor()Default constructor using Van Leer limiter. -
MUSCLReconstructor
Constructor with specified limiter.- Parameters:
limiter- Slope limiter type
-
-
Method Details
-
reconstruct
public MUSCLReconstructor.ReconstructedPair reconstruct(double uLeft, double uCenter, double uRight, double uFarRight) Reconstruct interface values between cells i and i+1.- Parameters:
uLeft- Value at cell i-1 (or i for boundary)uCenter- Value at cell iuRight- Value at cell i+1uFarRight- Value at cell i+2 (or i+1 for boundary)- Returns:
- Reconstructed left and right states at interface i+1/2
-
reconstructArray
Reconstruct interface values for an array of cell values.Returns arrays of left and right states at each interface.
- Parameters:
u- Cell-average values [0..n-1]- Returns:
- Array of reconstructed pairs at interfaces [0..n-2]
-
calcLimitedSlope
public double calcLimitedSlope(double deltaLeft, double deltaRight) Calculate limited slope using current limiter.- Parameters:
deltaLeft- Left difference (U_i - U_{i-1})deltaRight- Right difference (U_{i+1} - U_i)- Returns:
- Limited slope
-
calcLimiter
public double calcLimiter(double r) Calculate slope limiter value φ(r).- Parameters:
r- Gradient ratio- Returns:
- Limiter value φ
-
minmod
public double minmod(double r) Minmod limiter: φ(r) = max(0, min(1, r)).Most diffusive limiter, very robust for shocks.
- Parameters:
r- Gradient ratio- Returns:
- Limiter value
-
vanLeer
public double vanLeer(double r) Van Leer limiter: φ(r) = (r + |r|) / (1 + |r|).Smooth limiter with good balance of accuracy and stability.
- Parameters:
r- Gradient ratio- Returns:
- Limiter value
-
vanAlbada
public double vanAlbada(double r) Van Albada limiter: φ(r) = (r² + r) / (r² + 1).C1 continuous, good for smooth solutions.
- Parameters:
r- Gradient ratio- Returns:
- Limiter value
-
superbee
public double superbee(double r) Superbee limiter: φ(r) = max(0, min(2r, 1), min(r, 2)).Least diffusive, compresses discontinuities but can overshoot.
- Parameters:
r- Gradient ratio- Returns:
- Limiter value
-
mc
public double mc(double r) MC (Monotonized Central) limiter: φ(r) = max(0, min(2, 2r, (1+r)/2)).Between minmod and superbee in diffusivity.
- Parameters:
r- Gradient ratio- Returns:
- Limiter value
-
minmod3
public double minmod3(double a, double b, double c) Classic minmod function for three arguments.- Parameters:
a- First valueb- Second valuec- Third value- Returns:
- Minmod result
-
getLimiterType
Get current limiter type.- Returns:
- Current slope limiter
-
setLimiterType
Set limiter type.- Parameters:
limiterType- New slope limiter
-
isSecondOrder
public boolean isSecondOrder()Check if reconstruction is second-order.- Returns:
- true if using a limiter (second-order)
-