Ejemplo n.º 1
0
def test_exact_neumann_dirichlet():
    L = 10

    N = 128
    Ns = 200 + 1 #20000 + 1
    algo = 1
    scheme = 1

    W, u0, x = init_chebyshev_fredrikson(N, L)
    u0[0] = 0.

    print 'ETDRK4 N = ', N, ' Ns = ', Ns-1
    #q3, x3 = cheb_mde_neumann_dirichlet_etdrk4(W, u0, L, Ns, algo, scheme)
    lbc = BC('Neumann')
    rbc = BC('Dirichlet')
    etdrk4_solver = ETDRK4(L, N, Ns, h=None, lbc=lbc, rbc=rbc)
    q3, x3 = etdrk4_solver.solve(W, u0)
    Q3 = 0.5 * L * cheb_quadrature_clencurt(q3)
    if scheme == 0:
        #data_name = 'benchmark/NBC-DBC/exact/ETDRK4_Cox_N' 
        data_name = 'ETDRK4_Cox_N' 
        data_name = data_name + str(N) + '_Ns' + str(Ns-1)
    else:
        #data_name = 'benchmark/NBC-DBC/exact/ETDRK4_Krogstad_N' 
        data_name = 'ETDRK4_Krogstad_N' 
        data_name = data_name + str(N) + '_Ns' + str(Ns-1)
    savemat(data_name,{
            'N':N, 'Ns':Ns-1, 'W':W, 'u0':u0, 'Lz':L,
            'x':x, 'q':q3, 'Q':Q3})
    plt.plot(x3, q3, 'r')
    plt.axis([0, 10, 0, 1.15])
    plt.xlabel('$z$')
    plt.ylabel('$q(z)$')
    plt.savefig(data_name, bbox_inches='tight')
    plt.show()
Ejemplo n.º 2
0
    def build_grid(self):
        config = self.config
        Lx = config.grid.Lx
        self.Lx = Lx
        La = config.uc.a
        self.La = La
        Ms = config.grid.Ms

        if os.path.exists(config.grid.field_data):
            print "Reading from ", config.grid.field_data
            mat = loadmat(config.grid.field_data)
            self.w = np.zeros(Lx, dtype=np.complex128)
            self.w = mat['w'].reshape(Lx)
        else:
            self.w = np.zeros(Lx, dtype=np.complex128)
            self.w.real = np.random.rand(Lx)
            self.w.imag = np.random.rand(Lx)
        self.phi = np.zeros(Lx, dtype=np.complex128)

        self.q = np.zeros((Ms, Lx), dtype=np.complex128)
        self.q[0].real = 1.

        ds = 1. / (Ms - 1)
        self.ds = ds
        lbc = BC(self.lbc, self.lbc_vc)
        rbc = BC(self.rbc, self.rbc_vc)

        self.q_solver = ETDRK4(La, Lx - 1, Ms, h=ds, lbc=lbc, rbc=rbc)

        self.lam = config.grid.lam[0]
Ejemplo n.º 3
0
    def build_grid(self):
        config = self.config
        Lx = config.grid.Lx
        self.Lx = Lx
        La = config.uc.a
        self.La = La
        Ms = config.grid.Ms

        if os.path.exists(config.grid.field_data):
            mat = loadmat(config.grid.field_data)
            self.w = mat['w']
        else:
            #self.w = np.zeros([Lx])
            self.w = np.random.rand(Lx)

        self.q = np.zeros((Ms, Lx))
        self.q[0, :] = 1.

        ds = 1. / (Ms - 1)
        self.ds = ds
        lbc = BC(self.lbc, self.lbc_vc)
        rbc = BC(self.rbc, self.rbc_vc)

        self.q_solver = ETDRK4(La, Lx - 1, Ms, h=ds, lbc=lbc, rbc=rbc)

        self.lam = config.grid.lam[0]
Ejemplo n.º 4
0
    def build_grid(self):
        config = self.config
        Lx = config.grid.Lx
        self.Lx = Lx
        La = config.uc.a
        self.La = La
        Ms = config.grid.Ms
        MsA = config.grid.vMs[0]
        MsB = config.grid.vMs[1]
        MsC = config.grid.vMs[2]

        if os.path.exists(config.grid.field_data):
            mat = loadmat(config.grid.field_data)
            self.wA = mat['wA']
            self.wB = mat['wB']
            self.wC = mat['wC']
        else:
            self.wA = np.random.rand(Lx)
            self.wB = np.random.rand(Lx)
            self.wC = np.random.rand(Lx)

        self.qA = np.zeros((MsA, Lx))
        self.qA[0, :] = 1.
        self.qAc = np.zeros((MsA, Lx))
        self.qB = np.zeros((MsB, Lx))
        self.qBc = np.zeros((MsB, Lx))
        self.qBc[0, :] = 1.
        self.qC = np.zeros((MsC, Lx))
        self.qC[0, :] = 1.
        self.qCc = np.zeros((MsC, Lx))

        ds = 1. / (Ms - 1)
        self.ds = ds
        lbcA = BC(self.lbc, self.lbc_vc)
        rbcA = BC(self.rbc, self.rbc_vc)
        self.kaA = lbcA.beta
        self.kbA = rbcA.beta
        fA = self.fA
        fB = self.fB
        lbcB = copy.deepcopy(lbcA)
        lbcB.beta = -lbcA.beta * fA / fB  # See Fredrickson Book p.194
        rbcB = copy.deepcopy(rbcA)
        rbcB.beta = -rbcA.beta * fA / fB
        lbcC = copy.deepcopy(lbcA)
        rbcC = copy.deepcopy(rbcA)

        self.qA_solver = ETDRK4(La, Lx - 1, MsA, h=ds, lbc=lbcA, rbc=rbcA)
        self.qAc_solver = ETDRK4(La, Lx - 1, MsA, h=ds, lbc=lbcA, rbc=rbcA)
        self.qB_solver = ETDRK4(La, Lx - 1, MsB, h=ds, lbc=lbcB, rbc=rbcB)
        self.qBc_solver = ETDRK4(La, Lx - 1, MsB, h=ds, lbc=lbcB, rbc=rbcB)
        self.qC_solver = ETDRK4(La, Lx - 1, MsC, h=ds, lbc=lbcC, rbc=rbcC)
        #self.qCc_solver = ETDRK4(La, Lx-1, MsC-1, h=ds, lbc=lbcC, rbc=rbcC)
        self.qCc_solver = ETDRK4(La, Lx - 1, MsC, h=ds, lbc=lbcC, rbc=rbcC)

        self.lamA = config.grid.lam[0]
        self.lamB = config.grid.lam[1]
        self.lamC = config.grid.lam[2]
        self.lamY = config.grid.lam[3]
        #self.yita = self.lamY  # for compressible model
        self.yita = np.zeros(Lx)  # for incompressible model
Ejemplo n.º 5
0
 def build_grid(self):
     config = self.config
     Nt = config.grid.Lx
     Nr = config.grid.Ly # Nr is the number of grid points, EVEN only.
     Nr2 = Nr / 2
     self.Nt = Nt
     self.Nr = Nr
     self.Nr2 = Nr2
     R = config.uc.a
     self.R = R
     Ms = config.grid.Ms
     MsA = config.grid.vMs[0]
     MsB = config.grid.vMs[1]
     # Initiate fields from file or random numbers
     if os.path.exists(config.grid.field_data):
         mat = loadmat(config.grid.field_data)
         self.wA = mat['wA']
         self.wB = mat['wB']
     else:
         self.wA = np.random.rand(Nt,Nr2) - 0.5
         self.wB = np.random.rand(Nt,Nr2) - 0.5
     self.qA = np.zeros((MsA, Nt, Nr2))
     self.qA[0,:,:] = 1.
     self.qAc = np.zeros((MsA, Nt, Nr2))
     self.qB = np.zeros((MsB, Nt, Nr2))
     self.qBc = np.zeros((MsB, Nt, Nr2))
     self.qBc[0,:,:] = 1.
     ds = 1. / (Ms - 1)
     self.ds = ds
     lbcA = BC(self.lbc, self.lbc_vc)
     rbcA = BC(self.rbc, self.rbc_vc)
     self.kaA = lbcA.beta; self.kbA = rbcA.beta
     fA = self.fA; fB = 1 - fA
     lbcB = copy.deepcopy(lbcA)
     lbcB.beta = -lbcA.beta * fA / fB # See Fredrickson Book p.194
     rbcB = copy.deepcopy(rbcA)
     rbcB.beta = -rbcA.beta * fA / fB
     print 'lbcA=', lbcA.__dict__
     print 'rbcA=', rbcA.__dict__
     print 'lbcB=', lbcB.__dict__
     print 'rbcB=', rbcB.__dict__
     self.qA_solver = ETDRK4Polar(R, Nr-1, Nt,
                                          MsA, h=ds, lbc=lbcA, rbc=rbcA)
     self.qAc_solver = ETDRK4Polar(R, Nr-1, Nt,
                                           MsA, h=ds, lbc=lbcA, rbc=rbcA)
     self.qB_solver = ETDRK4Polar(R, Nr-1, Nt,
                                          MsB, h=ds, lbc=lbcB, rbc=rbcB)
     self.qBc_solver = ETDRK4Polar(R, Nr-1, Nt,
                                           MsB, h=ds, lbc=lbcB, rbc=rbcB)
     self.lamA = config.grid.lam[0]
     self.lamB = config.grid.lam[1]
     self.lamY = config.grid.lam[2]
     self.yita = self.lamY # For compressible model
Ejemplo n.º 6
0
 def build_grid(self):
     config = self.config
     Lx = config.grid.Lx
     Ly = config.grid.Ly
     self.Lx = Lx
     self.Ly = Ly
     La = config.uc.a
     Lb = config.uc.b
     self.La = La
     self.Lb = Lb
     Ms = config.grid.Ms
     MsA = config.grid.vMs[0]
     MsB = config.grid.vMs[1]
     if os.path.exists(config.grid.field_data):
         mat = loadmat(config.grid.field_data)
         self.wA = mat['wA']
         self.wB = mat['wB']
     else:
         self.wA = np.random.rand(Lx,Ly) - 0.5
         self.wB = np.random.rand(Lx,Ly) - 0.5
     self.qA = np.zeros((MsA, Lx, Ly))
     self.qA[0,:,:] = 1.
     self.qAc = np.zeros((MsA, Lx, Ly))
     self.qB = np.zeros((MsB, Lx, Ly))
     self.qBc = np.zeros((MsB, Lx, Ly))
     self.qBc[0,:,:] = 1.
     ds = 1. / (Ms - 1)
     self.ds = ds
     lbcA = BC(self.lbc, self.lbc_vc)
     rbcA = BC(self.rbc, self.rbc_vc)
     self.kaA = lbcA.beta; self.kbA = rbcA.beta
     fA = self.fA; fB = 1 - fA
     lbcB = copy.deepcopy(lbcA)
     lbcB.beta = -lbcA.beta * fA / fB # See Fredrickson Book p.194
     rbcB = copy.deepcopy(rbcA)
     rbcB.beta = -rbcA.beta * fA / fB
     print 'lbcA=', lbcA.__dict__
     print 'rbcA=', rbcA.__dict__
     print 'lbcB=', lbcB.__dict__
     print 'rbcB=', rbcB.__dict__
     self.qA_solver = ETDRK4FxCy(La, Lb, Lx, Ly-1,
                                          MsA, h=ds, lbc=lbcA, rbc=rbcA)
     self.qAc_solver = ETDRK4FxCy(La, Lb, Lx, Ly-1,
                                           MsA, h=ds, lbc=lbcA, rbc=rbcA)
     self.qB_solver = ETDRK4FxCy(La, Lb, Lx, Ly-1,
                                          MsB, h=ds, lbc=lbcB, rbc=rbcB)
     self.qBc_solver = ETDRK4FxCy(La, Lb, Lx, Ly-1,
                                           MsB, h=ds, lbc=lbcB, rbc=rbcB)
     self.lamA = config.grid.lam[0]
     self.lamB = config.grid.lam[1]
     self.lamY = config.grid.lam[2]
     self.yita = self.lamY # For compressible model
Ejemplo n.º 7
0
    def __init__(self,
                 Lx,
                 Ly,
                 Nx,
                 Ny,
                 Ns,
                 h=None,
                 c=1.0,
                 lbc=BC(),
                 rbc=BC(),
                 algo=1,
                 scheme=1):
        '''
        The PDE is in 2D,
            du/dt = cLu - wu
        where u=u(x,y), L=d^2/dx^2 + d^2/dy^2, w=w(x,y), c is a constant.
        First, perform a FFT in x direction to obtain
            du(kx,y)/dt = c L u(kx,y) - {kx^2 u(kx,y) + Fx[w(x,y)u(x,y)]}
        where L = D^2, with D^2 the Chebyshev 2nd order differential matrix,
        and kx^2 the d^2/dx^2 in Fourier space, see detail in
        the Notebook (page 2014.5.5).

        The defaut left BC and right BC are DBCs.

        Test: PASSED, 2014.5.8. (However, very small ds is required!)

        :param:Lx: physical size of the 1D spacial grid.
        :param:Lx: physical size of the 1D spacial grid.
        :param:Ns: number of grid points in time.
        :param:lbc: left boundary condition.
        :type:lbc: class BC
        :param:rbc: right boundary condition.
        :type:rbc: class BC
        :param:h: time step.
        :param:algo: algorithm for calculation of RK4 coefficients.
        :param:scheme: RK4 scheme.
        '''
        self.Lx, self.Ly = Lx, Ly
        self.Nx, self.Ny = Nx, Ny
        self.Ns = Ns
        if h is None:
            self.h = 1. / (Ns - 1)
        else:
            self.h = h
        self.c = c
        self.lbc = lbc
        self.rbc = rbc
        self.algo = algo
        self.scheme = scheme

        self.update()
Ejemplo n.º 8
0
 def build_grid(self):
     config = self.config
     Lx = config.grid.Lx
     self.Lx = Lx
     L = config.uc.a  # in unit of Rg
     self.L = L / self.z_hat  # in unit of \hat{z}
     N = Lx - 1
     self.N = N
     Ms = config.grid.Ms
     if os.path.exists(config.grid.field_data):
         mat = loadmat(config.grid.field_data)
         self.w = mat['w']
         self.w.shape = (self.w.size, )
     else:
         self.w = np.random.rand(Lx)
     self.q = np.zeros((Ms, Lx))
     self.q[0, :] = 1.
     self.qc = np.zeros((Ms, Lx))
     lbc = BC(self.lbc, self.lbc_vc)
     rbc = BC(self.rbc, self.rbc_vc)
     h = 1. / (Ms - 1)
     c = 0.25 / self.beta
     self.q_solver = ETDRK4(L,
                            N,
                            Ms,
                            h=h,
                            c=c,
                            lbc=lbc,
                            rbc=rbc,
                            algo=1,
                            scheme=1)
     #self.qc_solver = ETDRK4(L, N, Ms, h=h, c=c, lbc=lbc, rbc=rbc,
     #                        algo=1, scheme=1)
     self.qc_solver = ETDRK4(L,
                             N,
                             Ms - 1,
                             h=h,
                             c=c,
                             lbc=lbc,
                             rbc=rbc,
                             algo=1,
                             scheme=1)  # CKE
Ejemplo n.º 9
0
def test_etdrk4_complex():
    '''
    The test case is according to R. C. Daileda Lecture notes.
            du/dt = (1/25) u_xx , x@(0,3)
    with boundary conditions:
            u(0,t) = 0
            u_x(3,t) = -(1/2) u(3,t)
            u(x,0) = 100*(1-x/3)
    Conclusion:
        We find that the numerical solution is much more accurate than the five
        term approximation of the exact analytical solution.
    '''
    Nx = 64
    Lx = 3
    t = 1.
    Ns = 101
    ds = t / (Ns - 1)

    ii = np.arange(Nx + 1)
    x = np.cos(np.pi * ii / Nx)  # yy [-1, 1]
    x = 0.5 * (x + 1) * Lx  # mapping to [0, Ly]
    w = np.zeros(Nx + 1, dtype=np.complex128)
    q = np.zeros([Ns, Nx + 1], dtype=np.complex128)
    q[0] = 100 * (1 - x / 3)
    # The approximation of exact solution by first 5 terms
    q_exact = 47.0449*np.exp(-0.0210*t)*np.sin(0.7249*x) + \
              45.1413*np.exp(-0.1113*t)*np.sin(1.6679*x) + \
              21.3586*np.exp(-0.2872*t)*np.sin(2.6795*x) + \
              19.3403*np.exp(-0.5505*t)*np.sin(3.7098*x) + \
              12.9674*np.exp(-0.9015*t)*np.sin(4.7474*x)
    lbc = BC(DIRICHLET, [0, 1, 0])
    rbc = BC(ROBIN, [1., 0.5, 0])

    q_solver = ETDRK4(Lx, Nx, Ns, h=ds, c=1. / 25, lbc=lbc, rbc=rbc)
    q1, x = q_solver.solve(w, q[0], q)
    plt.plot(x, q[0].real, label='q_0')
    plt.plot(x, q1.real, label='q_solution1')
    plt.plot(x, q[-1].real, label='q_solution2')
    plt.plot(x, q_exact, label='q_exact')
    plt.legend(loc='best')
    plt.show()
Ejemplo n.º 10
0
    def build_grid(self):
        config = self.config
        Lx = config.grid.Lx
        self.Lx = Lx
        La = config.uc.a
        self.La = La
        Ms = config.grid.Ms

        if os.path.exists(config.grid.field_data):
            print "Reading from ", config.grid.field_data
            mat = loadmat(config.grid.field_data)
            self.w = np.zeros(Lx, dtype=np.complex128)
            #try:
            #    self.w = mat['iw_avg'][0,:]
            #except:
            #    self.w.real = mat['w'].reshape(Lx)
            #y = np.arange(-1, 1, 2.0/Lx)
            #self.w.real = cheb_interpolation_1d(y, w.real)
            #self.w.imag = cheb_interpolation_1d(y, w.imag)
            self.phi = np.zeros(Lx, dtype=np.complex128)
            try:
                self.phi.real = mat['phi_avg'][0, :].real
            except:
                self.phi.real = mat['phi'][0, :]
            #self.phi.real = cheb_interpolation_1d(y, phi.real)
        else:
            print "Abort: you must provide a density field."
            exit(1)

        self.q = np.zeros((Ms, Lx), dtype=np.complex128)
        self.q[0].real = 1.

        ds = 1. / (Ms - 1)
        self.ds = ds
        lbc = BC(self.lbc, self.lbc_vc)
        rbc = BC(self.rbc, self.rbc_vc)

        self.q_solver = ETDRK4(La, Lx - 1, Ms, h=ds, lbc=lbc, rbc=rbc)

        self.lam = config.grid.lam[0]
Ejemplo n.º 11
0
    def __init__(self, Lx, N, Ns, h=None, bc=BC()):
        '''
        Both BCs are DBCs or NBCs. The default is DBC.

        :param:Lx: physical size of the 1D spacial grid.
        :param:Ns: number of grid points in time.
        :param:N: number of grid points in space.
        :param:h: time step.
        '''
        self.Lx = Lx
        self.N = N
        self.Ns = Ns
        self.bc = bc
        if h is None:
            self.h = 1. / (Ns - 1)
        else:
            self.h = h

        self.update()
Ejemplo n.º 12
0
def test_etdrk4polar():
    '''
    The PDE is
            du/dt = (d^2/dr^2 + (1/r)d/dr + (1/r^2)d^2/dtheta^2) u - w u
    in the domain [0,R]x[0,2*pi] for time t=0 to t=1,
    with boundary conditions
            u(r,theta,t) = u(r,theta+2*pi,t) # periodic in theta direction
            du/dr |{r=1} + ka u(r=1) = 0, # RBC
    '''
    Nt = 64  # theta
    R = 1.0  # r [0, R]
    Nr = 63  # Nr must be odd
    N2 = (Nr + 1) / 2
    Ns = 101
    ds = 1. / (Ns - 1)

    # Periodic in x direction, Fourier
    tt = np.arange(Nt + 1) * 2 * np.pi / Nt
    # Non-periodic in y direction, Chebyshev
    ii = np.arange(Nr + 1)
    rr = np.cos(np.pi * ii / Nr)  # rr [-1, 1]
    rr = rr[:N2]  # rr in (0, 1] with r[0] = 1
    w = np.zeros([Nt, N2])
    A = 1.0
    q = np.zeros([Ns, Nt, N2])
    r, t = np.meshgrid(rr, tt[:-1])
    rp, tp = np.meshgrid(rr, tt)
    q[0] = np.exp(-A * r**2) * np.sin(t)
    w = np.cos(t)
    #plt.contourf(r*np.cos(t), r*np.sin(t), q[0])
    q0p = np.zeros([Nt + 1, N2])
    q0p[:-1, :] = q[0]
    q0p[-1, :] = q[0, 0, :]
    print rp.shape, tp.shape, q0p.shape
    scft_contourf(rp * np.cos(tp), rp * np.sin(tp), q0p)
    plt.xlabel('q[0]')
    plt.show()
    #plt.contourf(r*np.cos(t), r*np.sin(t), w)
    wp = np.zeros([Nt + 1, N2])
    wp[:-1, :] = w
    wp[-1, :] = w[0, :]
    scft_contourf(rp * np.cos(tp), rp * np.sin(tp), wp)
    plt.xlabel('w')
    plt.show()

    #lbc = BC(DIRICHLET, [0.0, 1.0, 0.0])
    #rbc = BC(DIRICHLET, [0.0, 1.0, 0.0])
    lbc = BC(ROBIN, [1.0, -0.1, 0.0])
    rbc = BC(ROBIN, [1.0, 0.1, 0.0])
    q_solver = ETDRK4Polar(R, Nr, Nt, Ns, h=ds, lbc=lbc, rbc=rbc)
    q1 = q_solver.solve(w, q[0], q)

    #plt.contourf(r*np.cos(t), r*np.sin(t), q1)
    q1p = np.zeros([Nt + 1, N2])
    q1p[:-1, :] = q1
    q1p[-1, :] = q1[0, :]
    scft_contourf(rp * np.cos(tp), rp * np.sin(tp), q1p)
    plt.xlabel('q_solution')
    plt.show()
    plt.plot(rr, q1[0, :], label='$\\theta=0$')
    plt.plot(rr, q1[Nt / 8, :], label='$\\theta=\pi/4$')
    plt.plot(rr, q1[Nt / 4, :], label='$\\theta=\pi/2$')
    plt.plot(rr, q1[3 * Nt / 8, :], label='$\\theta=3\pi/4$')
    plt.plot(rr, q1[Nt / 2, :], label='$\\theta=\pi$')
    plt.plot(rr, q1[5 * Nt / 8, :], label='$\\theta=5\pi/4$')
    plt.plot(rr, q1[3 * Nt / 4, :], label='$\\theta=3\pi/2$')
    plt.plot(rr, q1[7 * Nt / 8, :], label='$\\theta=7\pi/4$')
    plt.legend(loc='best')
    plt.show()
Ejemplo n.º 13
0
def test_exact_neumann(osc=0,oscheb=0,etdrk4=0):
    L = 10.0

    if osc:
        N = 128
        Ns = 1000 + 1 #20000 + 1
        W, u0, x = init_fourier(N, L)
        print 'OSC N = ', N, ' Ns = ', Ns-1
        #q1, x1 = cheb_mde_osc(W, u0, L, Ns)
        osc_solver = OSC(L, N, Ns)
        q1, x1 = osc_solver.solve(W, u0)
        Q1 = L * simps(q1, dx=1./N)
        #data_name = 'benchmark/NBC-NBC/exact/OSS_N' 
        data_name = 'OSS_N' 
        data_name = data_name + str(N) + '_Ns' + str(Ns-1)
        savemat(data_name,{
                'N':N, 'Ns':Ns-1, 'W':W, 'u0':u0, 'Lz':L,
                'x':x, 'q':q1, 'Q':Q1})
        plt.plot(x1, q1, 'b')
        plt.axis([0, 10, 0, 1.15])
        plt.xlabel('$z$')
        plt.ylabel('$q(z)$')
        plt.savefig(data_name, bbox_inches='tight')
        #plt.show()

    if oscheb:
        N = 128
        Ns = 200 + 1 #20000 + 1
        W, u0, x = init_chebyshev_fredrikson(N, L)
        print 'OSCHEB N = ', N, ' Ns = ', Ns-1
        #q2 = cheb_mde_neumann_oscheb(W, u0, L, Ns)
        oscheb_sovler = OSCHEB(L, N, Ns, bc=BC('Neumann'))
        q2, x2 = oscheb_sovler.solve(W, u0)
        Q2 = 0.5 * L * cheb_quadrature_clencurt(q2)
        #data_name = 'benchmark/NBC-NBC/exact/OSCHEB_N'
        data_name = 'OSCHEB_N'
        data_name = data_name + str(N) + '_Ns' + str(Ns-1)
        savemat(data_name,{
                'N':N, 'Ns':Ns-1, 'W':W, 'u0':u0, 'Lz':L,
                'x':x2, 'q':q2, 'Q':Q2})
        plt.plot(x2, q2, 'g')
        plt.axis([0, 10, 0, 1.15])
        plt.xlabel('$z$')
        plt.ylabel('$q(z)$')
        plt.savefig(data_name, bbox_inches='tight')
        #plt.show()

    if etdrk4:
        N = 128
        Ns = 200 + 1
        algo = 1
        scheme = 1
        W, u0, x = init_chebyshev_fredrikson(N, L)
        print 'ETDRK4 N = ', N, ' Ns = ', Ns-1
        #q3, x3 = cheb_mde_neumann_etdrk4(W, u0, L, Ns, None, algo, scheme)
        lbc = BC('Neumann')
        rbc = BC('Neumann')
        etdrk4_solver = ETDRK4(L, N, Ns, h=None, lbc=lbc, rbc=rbc)
        q3, x3 = etdrk4_solver.solve(W, u0)
        Q3 = 0.5 * L * cheb_quadrature_clencurt(q3)
        #if scheme == 0:
        #    data_name = 'benchmark/NBC-NBC/exact/ETDRK4_Cox_N' 
        #    data_name = data_name + str(N) + '_Ns' + str(Ns-1)
        #else:
        #    data_name = 'benchmark/NBC-NBC/exact/ETDRK4_Krogstad_N' 
        #    data_name = data_name + str(N) + '_Ns' + str(Ns-1)
        #savemat(data_name,{
        #        'N':N, 'Ns':Ns-1, 'W':W, 'u0':u0, 'Lz':L,
        #        'x':x, 'q':q3, 'Q':Q3})
        plt.plot(x3, q3, 'r')
        plt.axis([0, 10, 0, 1.15])
        plt.xlabel('$z$')
        plt.ylabel('$q(z)$')
        #plt.savefig(data_name, bbox_inches='tight')
        plt.show()
Ejemplo n.º 14
0
def test_etdrk4fxcy():
    '''
    The test function is
            u = e^[f(x,y) - t]
    where
            f(x,y) = h(x) + g(y)
    Assume it is the solution of following PDE
            du/dt = (d^2/dx^2 + d^2/dy^2) u - w(x,y)u
    in the domain [0,Lx]x[0,Ly] for time t=0 to t=1,
    with boundary conditions
            u(x+Lx,y,t) = u(x,y,t) # periodic in x direction
            d/dy[u(x,y=0,t)] = ka u(y=0)
            d/dy[u(x,y=Ly,t)] = kb u(y=Ly)
    To generate a suitable solution, we assume
            h(x) = sin(x)
            h_x = dh/dx = cos(x)
            h_xx = d^2h/dx^2 = -sin(x)
    since it is periodic in x direction.
    The corresponding w(x,y) is
            w(x,y) = h_xx + g_yy + (h_x)^2 + (g_y)^2 + 1

    1. For homogeneous NBC (ka=kb=0), a suitable g(y) is
            g(y) = Ay^2(2y-3)/6
            g_y = A(y^2-y)  # g_y(y=0)=0, g_y(y=1)=0
            g_yy = A(2*y-1)
    where A is a positive constant.
        Lx = 2*pi, Ly = 1.0, Nx = 64, Ny =32, Ns = 21
    is a good parameter set. Note the time step ds = 1/(Ns-1) = 0.05 is very
    large.
    2. For homogeneous DBC, an approximate g(y) is
            g(y) = -A(y-1)^2
            g_y = -2A(y-1)
            g_yy = -2A
    where A is a positive and large constant.
        Lx = 2*pi, Ly = 2.0, Nx = 64, Ny =32, Ns = 101
    is a good parameter set.
    3. For RBC, g(y) is given by
            g(y) = -Ay
            g_y = -A # ka=kb=-A
            g_yy = 0
       A is a positive constant.
       Numerical result is different than the analytical one.
    '''
    Lx = 2 * np.pi  # x [0, Lx]
    Nx = 64
    Ly = 1.0  # y [0, Ly]
    Ny = 127
    Ns = 101
    ds = 1. / (Ns - 1)

    # Periodic in x direction, Fourier
    xx = np.arange(Nx) * Lx / Nx
    # Non-periodic in y direction, Chebyshev
    ii = np.arange(Ny + 1)
    yy = np.cos(np.pi * ii / Ny)  # yy [-1, 1]
    yy = 0.5 * (yy + 1) * Ly  # mapping to [0, Ly]
    w = np.zeros([Nx, Ny + 1])
    A = 1.0
    q = np.zeros([Ns, Nx, Ny + 1])
    q_exact = np.zeros([Nx, Ny + 1])
    #q[0] = 1.
    for i in xrange(Nx):
        for j in xrange(Ny + 1):
            x = xx[i]
            y = yy[j]
            # RBC
            #q_exact[i,j] = np.exp(-A*y + np.sin(x) - 1)
            #q[0,i,j] = np.exp(-A*y + np.sin(x))
            #w[i,j] = np.cos(x)**2 - np.sin(x) + A**2 + 1
            # homogeneous NBC
            q_exact[i, j] = np.exp(A * y**2 * (2 * y - 3) / 6 + np.sin(x) - 1)
            q[0, i, j] = np.exp(A * y**2 * (2 * y - 3) / 6 + np.sin(x))
            w[i, j] = (
                A * y *
                (y - 1))**2 + np.cos(x)**2 - np.sin(x) + A * (2 * y - 1) + 1
            # homogeneous DBC
            #q[0,i,j] = np.exp(-A*(y-1)**2 + np.sin(x))
            #q_exact[i,j] = np.exp(-A*(y-1)**2 + np.sin(x) + 1)
            #w[i, j] = np.cos(x)**2 - np.sin(x) + 4*A**2 + (2*A*(y-1))**2 + 1
            # Fredrickson
            #sech = 1. / np.cosh(0.25*(6*y[j]-3*Ly))
            #w[i,j] = (1 - 2*sech**2)*(np.sin(2*np.pi*x[i]/Lx)+1)
            #w[i,j] = (1 - 2*sech**2)

    x = xx
    y = yy
    plt.imshow(w)
    plt.xlabel('w')
    plt.show()
    plt.plot(x, w[:, Ny / 2])
    plt.xlabel('w(x)')
    plt.show()
    plt.plot(y, w[Nx / 4, :])
    plt.xlabel('w(y)')
    plt.show()

    # DBC
    #lbc = BC(DIRICHLET, [0.0, 1.0, 0.0])
    #rbc = BC(DIRICHLET, [0.0, 1.0, 0.0])
    # RBC
    #lbc = BC(ROBIN, [1.0, A, 0.0])
    #rbc = BC(ROBIN, [1.0, A, 0.0])
    # NBC
    lbc = BC(ROBIN, [1.0, 0, 0.0])
    rbc = BC(ROBIN, [1.0, 0, 0.0])
    #q_solver = ETDRK4FxCy(Lx, Ly, Nx, Ny, Ns, h=ds, lbc=lbc, rbc=rbc)
    q_solver = ETDRK4FxCy2(Lx, Ly, Nx, Ny, Ns, h=ds, lbc=lbc, rbc=rbc)
    M = 100  # Took 1117.6 x 4 seconds for cpu one core
    with Timer() as t:
        for m in xrange(M):
            q1 = q_solver.solve(w, q[0], q)
    print "100 runs took ", t.secs, " seconds."

    print 'Error =', np.max(np.abs(q1 - q_exact))

    plt.imshow(q[0])
    plt.xlabel('q_0')
    plt.show()
    plt.imshow(q1)
    plt.xlabel('q_solution')
    plt.show()
    plt.imshow(q_exact)
    plt.xlabel('q_exact')
    plt.show()
    plt.plot(x, q[0, :, Ny / 2], label='q0')
    plt.plot(x, q1[:, Ny / 2], label='q_solution')
    plt.plot(x, q_exact[:, Ny / 2], label='q_exact')
    plt.legend(loc='best')
    plt.xlabel('q[:,Ny/2]')
    plt.show()
    plt.plot(y, q[0, Nx / 4, :], label='q0')
    plt.plot(y, q1[Nx / 4, :], label='q_solution')
    plt.plot(y, q_exact[Nx / 4, :], label='q_exact')
    plt.legend(loc='best')
    plt.xlabel('q[Nx/4,:]')
    plt.show()
    plt.plot(y, q[0, Nx * 3 / 4, :], label='q0')
    plt.plot(y, q1[Nx * 3 / 4, :], label='q_solution')
    plt.plot(y, q_exact[Nx * 3 / 4, :], label='q_exact')
    plt.legend(loc='best')
    plt.xlabel('q[Nx*3/4,:]')
    plt.show()
    exit()

    # Check with ETDRK4
    sech = 1. / np.cosh(0.25 * (6 * y - 3 * Ly))
    w1 = 1 - 2 * sech**2
    plt.plot(y, w1)
    plt.show()
    q = np.zeros([Ns, Ny + 1])
    q[0] = 1.
    q_solver = ETDRK4(Ly, Ny, Ns, h=ds, lbc=lbc, rbc=rbc)
    q1, y = q_solver.solve(w1, q[0], q)
    plt.plot(y, q1)
    plt.show()
Ejemplo n.º 15
0
def test_etdrk4fxycz():
    '''
    Assume the following PDE in cylindrical coordinates,
            du/dt = (d^2/dx^2 + d^2/dy^2 + d^2/dz^2) u - w u
    with u = u(x,y,z), w=w(x,y,z) in the domain
            [0, Lx] x [0, Ly] x [0, Lz]
    for time t=0 to t=1,
    with boundary conditions
            u(x,y,z,t) = u(x+Lx,y,z,t) # periodic in x direction
            u(x,y,z,t) = u(x,y+Ly,z,t) # periodic in y direction
            d/dr[u(z=0,y,z,t)] = ka u(z=0)
            d/dr[u(z=Lz,y,z,t)] = kb u(z=Lz)
    Test: PASSED 2013.8.16.
    '''
    Lx = 2.0 * np.pi  # x [0, Lx]
    Nx = 32
    Ly = 2.0 * np.pi  # y [0, Ly]
    Ny = 32
    Lz = 3.0  # z [0, Lz]
    Nz = 63
    Np = 64
    Ns = 101
    ds = 1. / (Ns - 1)
    show3d = False

    # Periodic in x direction, Fourier
    xx = np.arange(Nx) * Lx / Nx
    yy = np.arange(Ny) * Ly / Ny
    # Non-periodic in z direction, Chebyshev
    ii = np.arange(Nz + 1)
    zz = np.cos(np.pi * ii / Nz)  # zz [-1, 1]
    zz = 0.5 * (zz + 1) * Lz  # mapping to [0, Lz]
    zzp = np.linspace(0, Lz, Np)
    x, y, z = np.meshgrid(xx, yy, zz)
    #xp, yp, zp = np.meshgrid(xx, yy, zzp)
    xp, yp, zp = np.mgrid[0:Lx * (1 - 1. / Nx):Nx * 1j,
                          0:Ly * (1 - 1. / Ny):Ny * 1j, 0:Lz:Np * 1j]
    A = 1.0
    w = np.sin(x + 2 * y) + np.cos(z)**2
    q = np.zeros([Ns, Nx, Ny, Nz + 1])
    q[0] = np.exp(np.sin(x * y) - A * z)
    q0p = np.exp(np.sin(xp * yp) - A * zp)  # the exact

    #q0 = q[0]
    #plt.imshow(q0[Nx/2,:,:])
    #plt.show()
    #plt.imshow(q0[:,Ny/2,:])
    #plt.show()
    #plt.imshow(q0[:,:,Nz/2])
    #plt.show()
    if show3d:
        q0p2 = cheb_interp3d_z(q[0], zzp)
        print np.mean(np.abs(q0p2 - q0p))  # the error of interpolation
        #plt.plot(zz, q[0, Nx/2, Ny/2], '.')
        #plt.plot(zzp, q0p[Nx/2, Ny/2])
        #plt.plot(zzp, q0p2[Nx/2, Ny/2], '^')
        #plt.show()
        mlab.clf()
        #s = mlab.pipeline.scalar_field(xp,yp,zp,q0p2)
        #mlab.pipeline.iso_surface(s,
        #contours=[q0p.min()+0.1*q0p.ptp(),],
        #                          opacity=0.3)
        mlab.contour3d(xp,
                       yp,
                       zp,
                       q0p2,
                       contours=16,
                       opacity=0.5,
                       transparent=True,
                       colormap='Spectral')
        mlab.show()

    #lbc = BC(DIRICHLET, [0.0, 1.0, 0.0])
    #rbc = BC(DIRICHLET, [0.0, 1.0, 0.0])
    lbc = BC(ROBIN, [1.0, 0.1, 0.0])
    rbc = BC(ROBIN, [1.0, -0.1, 0.0])
    q_solver = ETDRK4FxyCz(Lx, Ly, Lz, Nx, Ny, Nz, Ns, h=ds, lbc=lbc, rbc=rbc)
    print 'Initiate done!'
    M = 10
    # 10 runs took 131.901 seconds for cpu one core (32x32x32)
    with Timer() as t:
        for m in xrange(M):
            q1 = q_solver.solve(w, q[0], q)
    print "10 runs took", t.secs, " seconds."
    print 'Solve done!'

    #plt.imshow(q1[Nx/2,:,:])
    #plt.show()
    #plt.imshow(q1[:,Ny/2,:])
    #plt.show()
    #plt.imshow(q1[:,:,Nz/2])
    #plt.show()
    if show3d:
        q1p = cheb_interp3d_z(q1, zzp)
        mlab.clf()
        mlab.contour3d(xp, yp, zp, q1p, contours=16, transparent=True)
        #mlab.contour3d(q1)
        mlab.show()
Ejemplo n.º 16
0
def test_etdrk4cylind():
    '''
    Assume the following PDE in cylindrical coordinates,
            du/dt = (d^2/dr^2 + (1/r)d/dr + (1/r^2)d^2/d^theta + d^2/dz^2) u -
            w*u
    with u = u(r,theta,z), w=w(r,theta,z) in the domain 
            [0, R] x [0, 2pi] x [0, Lz]
    for time t=0 to t=1,
    with boundary conditions
            d/dr[u(r=R,theta,z,t)] = ka u(r=R)
            u(r,theta,z,t) = u(r,theta+2pi,z,t) # periodic in theta direction
            u(r,theta,z,t) = u(r,theta,z+lz,t) # periodic in z direction
    '''
    Lt = 2*np.pi
    Nt = 32 # theta
    Lz = 4.0
    Nz = 64
    R = 2.0 # r [0, R]
    Nr = 23 # Nr must be odd
    N2 = (Nr + 1) / 2
    Nxp = 20
    Nyp = 20
    Nzp = 15
    Ns = 101 
    ds = 1. / (Ns - 1)
    
    # Periodic in x and z direction, Fourier
    tt = np.arange(Nt) * Lt / Nt
    zz = np.arange(Nz) * Lz / Nz
    # Non-periodic in r direction, Chebyshev
    ii = np.arange(Nr+1)
    rr = np.cos(np.pi * ii / Nr) # rr [-1, 1]
    rr = rr[:N2] # rr in (0, 1] with r[0] = 1
    rrp = np.linspace(0, R, Nxp)
    A = 1.0
    q = np.zeros([Ns, Nt, Nz, N2])
    z, t, r = np.meshgrid(zz, tt, rr)
    zp, tp, rp = np.meshgrid(zz, tt, rrp)
    print tt.shape, zz.shape, rr.shape
    print t.shape, z.shape, r.shape
    print q[0].shape
    q[0] = r**2 * np.cos(t) + A*z
    xp = rp * np.cos(tp)
    yp = rp * np.sin(tp)
    print xp.max(), xp.min()
    print yp.max(), yp.min()
    print zp.max(), zp.min()
    q0 = cheb_interp3d_r(q[0], rrp)
    xpp, ypp, zpp = np.mgrid[-R:R:Nxp*1j,
                    -R:R:Nyp*1j,
                    0:Lz*(1-1./Nz):Nzp*1j]
    q0p = griddata((xp.ravel(),yp.ravel(),zp.ravel()), q0.ravel(), 
                   (xpp, ypp, zpp), method='linear')
    print q0p.shape
    mlab.contour3d(xpp, ypp, zpp, q0p, 
                   contours=64, transparent=True, colormap='Spectral')
    mlab.show()
    exit()

    q0p[:-1,:] = q[0]
    q0p[-1,:] = q[0,0,:]
    print rp.shape, tp.shape, q0p.shape
    scft_contourf(rp*np.cos(tp), rp*np.sin(tp), q0p)
    plt.xlabel('q[0]')
    plt.show()
    #plt.contourf(r*np.cos(t), r*np.sin(t), w)
    w = np.zeros([Nt, Nz, N2])
    w = np.cos(t)
    wp = np.zeros([Nt+1,N2])
    wp[:-1,:] = w
    wp[-1,:] = w[0,:]
    scft_contourf(rp*np.cos(tp), rp*np.sin(tp), wp)
    plt.xlabel('w')
    plt.show()

    #lbc = BC(DIRICHLET, [0.0, 1.0, 0.0]) 
    #rbc = BC(DIRICHLET, [0.0, 1.0, 0.0]) 
    lbc = BC(ROBIN, [1.0, -0.1, 0.0]) 
    rbc = BC(ROBIN, [1.0, 0.1, 0.0]) 
    q_solver = ETDRK4Cylind(R, Nr, Nt, Ns, h=ds, lbc=lbc, rbc=rbc)
    q1 = q_solver.solve(w, q[0], q)

    #plt.contourf(r*np.cos(t), r*np.sin(t), q1)
    q1p = np.zeros([Nt+1,N2])
    q1p[:-1,:] = q1
    q1p[-1,:] = q1[0,:]
    scft_contourf(rp*np.cos(tp), rp*np.sin(tp), q1p)
    plt.xlabel('q_solution')
    plt.show()
    plt.plot(rr, q1[0,:], label='$\\theta=0$')
    plt.plot(rr, q1[Nt/8,:], label='$\\theta=\pi/4$')
    plt.plot(rr, q1[Nt/4,:], label='$\\theta=\pi/2$')
    plt.plot(rr, q1[3*Nt/8,:], label='$\\theta=3\pi/4$')
    plt.plot(rr, q1[Nt/2,:], label='$\\theta=\pi$')
    plt.plot(rr, q1[5*Nt/8,:], label='$\\theta=5\pi/4$')
    plt.plot(rr, q1[3*Nt/4,:], label='$\\theta=3\pi/2$')
    plt.plot(rr, q1[7*Nt/8,:], label='$\\theta=7\pi/4$')
    plt.legend(loc='best')
    plt.show()