def compute_and_plot_solutions():
    """
    Routine to compute and plot the solutions of SDC(2), IMEX, BDF-2 and RK for a multiscale problem
    """

    num_procs = 1

    t0 = 0.0
    Tend = 3.0
    nsteps = 154  # 154 is value in Vater et al.
    dt = Tend / float(nsteps)

    # initialize level parameters
    level_params = dict()
    level_params['restol'] = 1E-10
    level_params['dt'] = dt

    # This comes as read-in for the step class
    step_params = dict()
    step_params['maxiter'] = 2

    # This comes as read-in for the problem class
    problem_params = dict()
    problem_params['cadv'] = 0.05
    problem_params['cs'] = 1.0
    problem_params['nvars'] = [(2, 512)]
    problem_params['order_adv'] = 5
    problem_params['waveno'] = 5

    # This comes as read-in for the sweeper class
    sweeper_params = dict()
    sweeper_params['collocation_class'] = CollGaussRadau_Right
    sweeper_params['num_nodes'] = 2

    # initialize controller parameters
    controller_params = dict()
    controller_params['logger_level'] = 30
    controller_params['hook_class'] = dump_energy

    # Fill description dictionary for easy hierarchy creation
    description = dict()
    description['problem_class'] = acoustic_1d_imex_multiscale
    description['problem_params'] = problem_params
    description['sweeper_class'] = imex_1st_order
    description['sweeper_params'] = sweeper_params
    description['step_params'] = step_params
    description['level_params'] = level_params

    # instantiate the controller
    controller = controller_nonMPI(num_procs=num_procs,
                                   controller_params=controller_params,
                                   description=description)

    # get initial values on finest level
    P = controller.MS[0].levels[0].prob
    uinit = P.u_exact(t0)

    # call main function to get things done...
    uend, stats = controller.run(u0=uinit, t0=t0, Tend=Tend)

    # instantiate standard integrators to be run for comparison
    trap = trapezoidal((P.A + P.Dx).astype('complex'), 0.5)
    bdf2_m = bdf2(P.A + P.Dx)
    dirk_m = dirk((P.A + P.Dx).astype('complex'), step_params['maxiter'])
    rkimex = rk_imex(P.A.astype('complex'), P.Dx.astype('complex'),
                     step_params['maxiter'])

    y0_tp = np.concatenate((uinit.values[0, :], uinit.values[1, :]))
    y0_bdf = y0_tp
    y0_dirk = y0_tp.astype('complex')
    y0_imex = y0_tp.astype('complex')

    # Perform time steps with standard integrators
    for i in range(0, nsteps):

        # trapezoidal rule step
        ynew_tp = trap.timestep(y0_tp, dt)

        # BDF-2 scheme
        if i == 0:
            ynew_bdf = bdf2_m.firsttimestep(y0_bdf, dt)
            ym1_bdf = y0_bdf
        else:
            ynew_bdf = bdf2_m.timestep(y0_bdf, ym1_bdf, dt)

        # DIRK scheme
        ynew_dirk = dirk_m.timestep(y0_dirk, dt)

        # IMEX scheme
        ynew_imex = rkimex.timestep(y0_imex, dt)

        y0_tp = ynew_tp
        ym1_bdf = y0_bdf
        y0_bdf = ynew_bdf
        y0_dirk = ynew_dirk
        y0_imex = ynew_imex

        # Finished running standard integrators
        unew_tp, pnew_tp = np.split(ynew_tp, 2)
        unew_bdf, pnew_bdf = np.split(ynew_bdf, 2)
        unew_dirk, pnew_dirk = np.split(ynew_dirk, 2)
        unew_imex, pnew_imex = np.split(ynew_imex, 2)

    fs = 8

    rcParams['figure.figsize'] = 2.5, 2.5
    # rcParams['pgf.rcfonts'] = False
    fig = plt.figure()

    sigma_0 = 0.1
    k = 7.0 * 2.0 * np.pi
    x_0 = 0.75
    x_1 = 0.25

    print('Maximum pressure in SDC: %5.3e' %
          np.linalg.norm(uend.values[1, :], np.inf))
    print('Maximum pressure in DIRK: %5.3e' %
          np.linalg.norm(pnew_dirk, np.inf))
    print('Maximum pressure in RK-IMEX: %5.3e' %
          np.linalg.norm(pnew_imex, np.inf))

    if dirk_m.order == 2:
        plt.plot(P.mesh,
                 pnew_bdf,
                 'd-',
                 color='c',
                 label='BDF-2',
                 markevery=(50, 75))
    p_slow = np.exp(
        -np.square(np.mod(P.mesh - problem_params['cadv'] * Tend, 1.0) - x_0) /
        (sigma_0 * sigma_0))
    plt.plot(P.mesh,
             p_slow,
             '--',
             color='k',
             markersize=fs - 2,
             label='Slow mode',
             dashes=(10, 2))
    if np.linalg.norm(pnew_imex, np.inf) <= 2:
        plt.plot(P.mesh,
                 pnew_imex,
                 '+-',
                 color='r',
                 label='IMEX(' + str(rkimex.order) + ')',
                 markevery=(1, 75),
                 mew=1.0)
    plt.plot(P.mesh,
             uend.values[1, :],
             'o-',
             color='b',
             label='SDC(' + str(step_params['maxiter']) + ')',
             markevery=(25, 75))
    plt.plot(P.mesh,
             pnew_dirk,
             '-',
             color='g',
             label='DIRK(' + str(dirk_m.order) + ')')

    plt.xlabel('x', fontsize=fs, labelpad=0)
    plt.ylabel('Pressure', fontsize=fs, labelpad=0)
    fig.gca().set_xlim([0, 1.0])
    fig.gca().set_ylim([-0.5, 1.1])
    fig.gca().tick_params(axis='both', labelsize=fs)
    plt.legend(loc='upper left',
               fontsize=fs,
               prop={'size': fs},
               handlelength=3)
    fig.gca().grid()
    filename = 'data/multiscale-K' + str(step_params['maxiter']) + '-M' + str(
        sweeper_params['num_nodes']) + '.png'
    plt.gcf().savefig(filename, bbox_inches='tight')
Пример #2
0
def compute_and_plot_dispersion():
    problem_params = dict()
    # SET VALUE FOR lambda_slow AND VALUES FOR lambda_fast ###
    problem_params['lambda_s'] = np.array([0.0])
    problem_params['lambda_f'] = np.array([0.0])
    problem_params['u0'] = 1.0

    # initialize sweeper parameters
    sweeper_params = dict()
    # SET TYPE AND NUMBER OF QUADRATURE NODES ###
    sweeper_params['collocation_class'] = CollGaussRadau_Right
    sweeper_params['do_coll_update'] = True
    sweeper_params['num_nodes'] = 3

    # initialize level parameters
    level_params = dict()
    level_params['dt'] = 1.0

    # fill description dictionary for easy step instantiation
    description = dict()
    description['problem_class'] = swfw_scalar  # pass problem class
    description['problem_params'] = problem_params  # pass problem parameters
    description['dtype_u'] = mesh  # pass data type for u
    description['dtype_f'] = rhs_imex_mesh  # pass data type for f
    description['sweeper_class'] = imex_1st_order  # pass sweeper
    description['sweeper_params'] = sweeper_params  # pass sweeper parameters
    description['level_params'] = level_params  # pass level parameters
    description['step_params'] = dict()  # pass step parameters

    # SET NUMBER OF ITERATIONS ###
    K = 3

    # ORDER OF DIRK/IMEX IS EQUAL TO NUMBER OF ITERATIONS AND THUS ORDER OF SDC ###
    dirk_order = K

    c_speed = 1.0
    U_speed = 0.05

    # now the description contains more or less everything we need to create a step
    S = step(description=description)

    L = S.levels[0]

    # u0 = S.levels[0].prob.u_exact(t0)
    # S.init_step(u0)
    QE = L.sweep.QE[1:, 1:]
    QI = L.sweep.QI[1:, 1:]
    Q = L.sweep.coll.Qmat[1:, 1:]
    nnodes = L.sweep.coll.num_nodes
    dt = L.params.dt

    Nsamples = 15
    k_vec = np.linspace(0, np.pi, Nsamples + 1, endpoint=False)
    k_vec = k_vec[1:]
    phase = np.zeros((3, Nsamples))
    amp_factor = np.zeros((3, Nsamples))

    for i in range(0, np.size(k_vec)):

        Cs = -1j * k_vec[i] * np.array([[0.0, c_speed], [c_speed, 0.0]],
                                       dtype='complex')
        Uadv = -1j * k_vec[i] * np.array([[U_speed, 0.0], [0.0, U_speed]],
                                         dtype='complex')

        LHS = np.eye(2 * nnodes) - dt * (np.kron(QI, Cs) + np.kron(QE, Uadv))
        RHS = dt * (np.kron(Q, Uadv + Cs) - np.kron(QI, Cs) -
                    np.kron(QE, Uadv))

        LHSinv = np.linalg.inv(LHS)
        Mat_sweep = np.linalg.matrix_power(LHSinv.dot(RHS), K)
        for k in range(0, K):
            Mat_sweep = Mat_sweep + np.linalg.matrix_power(LHSinv.dot(RHS),
                                                           k).dot(LHSinv)
        ##
        # ---> The update formula for this case need verification!!
        update = dt * np.kron(L.sweep.coll.weights, Uadv + Cs)

        y1 = np.array([1, 0], dtype='complex')
        y2 = np.array([0, 1], dtype='complex')
        e1 = np.kron(np.ones(nnodes), y1)
        stab_fh_1 = y1 + update.dot(Mat_sweep.dot(e1))
        e2 = np.kron(np.ones(nnodes), y2)
        stab_fh_2 = y2 + update.dot(Mat_sweep.dot(e2))
        stab_sdc = np.column_stack((stab_fh_1, stab_fh_2))

        # Stability function of backward Euler is 1/(1-z); system is y' = (Cs+Uadv)*y
        # stab_ie = np.linalg.inv( np.eye(2) - step.status.dt*(Cs+Uadv) )

        # For testing, insert exact stability function exp(-dt*i*k*(Cs+Uadv)
        # stab_fh = la.expm(Cs+Uadv)

        dirkts = dirk(Cs + Uadv, dirk_order)
        stab_fh1 = dirkts.timestep(y1, 1.0)
        stab_fh2 = dirkts.timestep(y2, 1.0)
        stab_dirk = np.column_stack((stab_fh1, stab_fh2))

        rkimex = rk_imex(M_fast=Cs, M_slow=Uadv, order=K)
        stab_fh1 = rkimex.timestep(y1, 1.0)
        stab_fh2 = rkimex.timestep(y2, 1.0)
        stab_rk_imex = np.column_stack((stab_fh1, stab_fh2))

        sol_sdc = findomega(stab_sdc)
        sol_dirk = findomega(stab_dirk)
        sol_rk_imex = findomega(stab_rk_imex)

        # Now solve for discrete phase
        phase[0, i] = sol_sdc.real / k_vec[i]
        amp_factor[0, i] = np.exp(sol_sdc.imag)
        phase[1, i] = sol_dirk.real / k_vec[i]
        amp_factor[1, i] = np.exp(sol_dirk.imag)
        phase[2, i] = sol_rk_imex.real / k_vec[i]
        amp_factor[2, i] = np.exp(sol_rk_imex.imag)

    rcParams['figure.figsize'] = 1.5, 1.5
    fs = 8
    fig = plt.figure()
    plt.plot(k_vec, (U_speed + c_speed) + np.zeros(np.size(k_vec)),
             '--',
             color='k',
             linewidth=1.5,
             label='Exact')
    plt.plot(k_vec,
             phase[1, :],
             '-',
             color='g',
             linewidth=1.5,
             label='DIRK(' + str(dirkts.order) + ')')
    plt.plot(k_vec,
             phase[2, :],
             '-+',
             color='r',
             linewidth=1.5,
             label='IMEX(' + str(rkimex.order) + ')',
             markevery=(2, 3),
             mew=1.0)
    plt.plot(k_vec,
             phase[0, :],
             '-o',
             color='b',
             linewidth=1.5,
             label='SDC(' + str(K) + ')',
             markevery=(1, 3),
             markersize=fs / 2)
    plt.xlabel('Wave number', fontsize=fs, labelpad=0.25)
    plt.ylabel('Phase speed', fontsize=fs, labelpad=0.5)
    plt.xlim([k_vec[0], k_vec[-1:]])
    plt.ylim([0.0, 1.1 * (U_speed + c_speed)])
    fig.gca().tick_params(axis='both', labelsize=fs)
    plt.legend(loc='lower left', fontsize=fs, prop={'size': fs - 2})
    plt.xticks([0, 1, 2, 3], fontsize=fs)
    filename = 'data/phase-K' + str(K) + '-M' + str(
        sweeper_params['num_nodes']) + '.png'
    plt.gcf().savefig(filename, bbox_inches='tight')

    fig = plt.figure()
    plt.plot(k_vec,
             1.0 + np.zeros(np.size(k_vec)),
             '--',
             color='k',
             linewidth=1.5,
             label='Exact')
    plt.plot(k_vec,
             amp_factor[1, :],
             '-',
             color='g',
             linewidth=1.5,
             label='DIRK(' + str(dirkts.order) + ')')
    plt.plot(k_vec,
             amp_factor[2, :],
             '-+',
             color='r',
             linewidth=1.5,
             label='IMEX(' + str(rkimex.order) + ')',
             markevery=(2, 3),
             mew=1.0)
    plt.plot(k_vec,
             amp_factor[0, :],
             '-o',
             color='b',
             linewidth=1.5,
             label='SDC(' + str(K) + ')',
             markevery=(1, 3),
             markersize=fs / 2)
    plt.xlabel('Wave number', fontsize=fs, labelpad=0.25)
    plt.ylabel('Amplification factor', fontsize=fs, labelpad=0.5)
    fig.gca().tick_params(axis='both', labelsize=fs)
    plt.xlim([k_vec[0], k_vec[-1:]])
    plt.ylim([k_vec[0], k_vec[-1:]])
    plt.legend(loc='lower left', fontsize=fs, prop={'size': fs - 2})
    plt.gca().set_ylim([0.0, 1.1])
    plt.xticks([0, 1, 2, 3], fontsize=fs)
    filename = 'data/ampfactor-K' + str(K) + '-M' + str(
        sweeper_params['num_nodes']) + '.png'
    plt.gcf().savefig(filename, bbox_inches='tight')