Carmen Guo
Explain why modelling advection is particularly difficult
Explain and contrast diffusion and dispersion
Define positive definite
[10]:
import context
from IPython.display import Image
import matplotlib.pyplot as plt
# make the plots happen inline
%matplotlib inline
# import the numpy array handling library
import numpy as np
# import the advection code from the numlabs directory
import numlabs.lab10.advection_funs as afs
The word advection means ‘transfer of heat through the horizontal motion of a flow’. More generally, we will consider a flow of some non-diffusive quantity. For example, consider the wind as a flow and the water vapour in the air as the non-diffusive quantity. Suppose the wind is travelling in the positive x direction, and we are considering the vapour concentration from
Assume that initially the distribution curve of the vapour is Gaussian (Figure Initial Distribution). Ideally, the water droplets move at the same speed as that of the air, so the distribution curve retains its initial shape as it travels along the x-axis. This process is described by the following PDE:
(Advection Eqn)
where
[11]:
Image(filename='images/initial.png',width='60%')
[11]:
Figure Initial Distribution: This is the initial distribution of water vapour concentration.
As you will see in the upcoming examples, it is not easy to obtain a satisfactory numerical solution to this PDE.
Let’s start off simple and solve the PDE (Advection Eqn) using centred differences, i.e., by expanding the time and spatial derivatives in the following way:
(Centered Difference Scheme)
where
The boundary conditions are :
The initial conditions are:
where
Now we need the values of
Substitution of the equations into the PDE yields
where
The function that computes the numerical solution using this scheme is in advection_funs.py. It is a python function advection(timesteps), which takes in the number of time steps as input and plots the distribution curve at 20 time steps.
We can see the problem with this scheme just by running the function with 10 time steps (Figure Distribution with Centered Scheme). Following is a plot of the distribution curve at the last time step.
[12]:
Image(filename='images/centered.png',width='60%')
[12]:
Figure Distribution with Centered Scheme: This is the distribution after 10 time steps approximated using the centred differencing scheme.
Comparing this curve with the initial state (Figure Initial Distribution), we can see that ripples are produced to the left of the curve which means negative values are generated. But water vapour does not have negative concentrations. The centred differencing scheme does not work well for the advection process because it is not positive definite, i.e., it generates negative values which are impossible in real life.
Another example of using the same scheme
[13]:
afs.advection(10, lab_example=True)
Let’s see what is wrong with our simple centred differencing scheme. We used centred differences to compute the time and spatial derivatives (Centered Difference Scheme). In other words,
If we use backward differences to approximate the spatial derivative but continue to use centred differences to approximate the time derivative, we will end up with an unstable scheme. Thus, we will use backward differences for both time and spatial derivatives. Now the time and spatial derivatives are given by:
(Upstream Scheme)
Substitution of the equations into the PDE yields
The boundary conditions and the initial conditions are the same as in the centred differencing scheme. This time we compute
The function that computes the solution using this scheme is in advection_funs.py. It is a python function advection2(timesteps), which takes in the number of time steps as input and plots the distribution curve at 20 time steps.
Although this scheme is positive definite and conservative (the area under the curve is the same as in the initial state), it introduces a new problem — diffusion. As the time step increases, you can see that the curve becomes wider and lower, i.e., it diffuses quickly (Figure Upstream Distribution).
[14]:
Image(filename='images/upstream.png',width='60%')
[14]:
Figure Upstream Distribution: This is the distribution after 60 time steps approximated using the upstream method.
But ideally, the curve should retain its original shape as the wave travels along the x-axis, so the upstream method is still not good enough for the advection problem. In the next section, we will present another method that does a better job.
Another example using the same scheme
[15]:
afs.advection2(60, lab_example=True)
In previous sections, we were concerned with values at grid points only, ie, values of
The PDE (Advection Eqn) is rewritten as :
(Flux Form Eqn)
where
Expand the time derivative using forward differences and the spatial derivative as follows:
where
Substituting the expanded derivatives into the PDE, we obtain the following recurrence formula (
Since flux is defined as the amount flowing through per unit time, we need to calculate the portion of the distribution curve in each grid box that passes the right boundary after
[16]:
Image(filename='images/flux.png',width='60%')
[16]:
Figure Amount Leaving Grid: After
As (Figure Amount Leaving Grid) shows, the distribution curve in grid box
Since we only know
where
The coefficients
Table :math:`ell = 0`:
Table :math:`ell = 1`: two representations
Table :math:`ell = 2`:
Table :math:`ell = 3`: two representations
Table :math:`ell = 4`:
Note that for even order polynomials, an odd number of
If we choose
Now we define
(Flux Leaving Eqn)
Note that we are integrating over
(Flux Eqn)
In this form, the scheme is conservative and weakly diffusive. But it still lacks positive definiteness. A sufficient condition for this is
(Positive Definite Condition)
That is, the total outflow is never negative and never greater than
[17]:
Image(filename='images/limit.png',width='60%')
[17]:
Figure Total in Cell: The shaded area is equal to
If the total outflow is larger than the shaded area
To satisfy the condition for positive definiteness (Positive Definite Condition) we need to guarantee that
(Normalization Eqn)
where
Since the total flow out of a grid box is always less than the total grid volume,
Now to satisfy the lower limit of the positive definiteness condition, (Positive Definite Condition) we need to make sure
If we are looking at the parts of the curve that are far away from the peak,
Combining all the conditions from above, the advection scheme is described as follows:
with
where
An example function for this scheme is in advection_funs.py. The python function advection3(timesteps, order) takes in 2 arguments, the first is the number of time steps to be computed, the second is the order of the polynomial for the approximation of
[18]:
afs.advection3(60, 4, lab_example=True)
[18]:
array([[0.00000000e+00, 0.00000000e+00, 7.44658307e-03, ...,
0.00000000e+00, 0.00000000e+00, 0.00000000e+00],
[0.00000000e+00, 0.00000000e+00, 0.00000000e+00, ...,
0.00000000e+00, 3.04613895e-05, 0.00000000e+00],
[0.00000000e+00, 0.00000000e+00, 0.00000000e+00, ...,
0.00000000e+00, 3.04613895e-05, 0.00000000e+00],
...,
[0.00000000e+00, 0.00000000e+00, 0.00000000e+00, ...,
0.00000000e+00, 3.04515565e-05, 9.83308063e-09],
[0.00000000e+00, 0.00000000e+00, 0.00000000e+00, ...,
0.00000000e+00, 3.04449199e-05, 1.64696507e-08],
[0.00000000e+00, 0.00000000e+00, 0.00000000e+00, ...,
0.00000000e+00, 3.04366818e-05, 2.47077779e-08]])
Using the Bott Scheme, modify initialize in advection_funs.py to solve the following advection problem: The wind is moving along the x-axis with speed u=20 m/s. The initial distribution curve is 50 km in width. Use your program to approximate the curve during 24 hours.
Run your program for different orders of approximating polynomials (up to 4). Compare the accuracy of approximation for different orders. Do you see better results with increasing order? Is this true for all orders from 0 to 4? Is there any particularity to odd and even order polynomials? What if you decrease the Courant number (value in front of dx/u for dt calculation)?
b) For odd ordered polynomials, advection_funs.py uses the representation of
The last scheme presented solves 2 problems introduced in the previous two schemes. The centred differencing scheme lacks positive definiteness because it is of second order accuracy, thus it introduces additional oscillations near the peak. The scheme presented here solves this problem by checking and normalising the relevant values (ie,
Experiments have shown that this scheme is numerically stable in most atmospheric situations. This scheme is only slightly unstable in the case of strong deformational flow field models.
For more detail about this advection scheme , please refer to (Bott, 1989). Since 1989, new advection schemes including the MUSCL and TVD have been developed and are more routinely used.
Bott, A., 1989: A positive definite advection scheme obtained by nonlinear renormalization of the advective fluxes. Monthly Weather Review, 117, 1006–1015.
[ ]: