예제 #1
0
def setup(test_functions, trial_functions, w_, w_1,
          dx, ds, normal,
          dirichlet_bcs, neumann_bcs,
          boundary_to_mark,
          permittivity,
          density, viscosity,
          solutes, enable_PF, enable_EC, enable_NS,
          surface_tension, dt, interface_thickness,
          grav_const, pf_mobility_coeff, pf_mobility,
          q_rhs,
          use_iterative_solvers,
          p_lagrange,
          **namespace):
    """ Set up problem. """
    # Constants
    sigma_bar = surface_tension*3./(2*math.sqrt(2))
    per_tau = df.Constant(1./dt)
    grav = df.Constant((0., -grav_const))
    gamma = pf_mobility_coeff
    eps = interface_thickness

    # Set up the fields
    funs_ = df.split(w_["NSPFEC"])
    funs_1 = df.split(w_1["NSPFEC"])
    field_number = 0
    if enable_NS:
        v = test_functions["NSPFEC"][field_number]
        u_ = funs_[field_number]
        u_1 = funs_1[field_number]
        field_number += 1

        q = test_functions["NSPFEC"][field_number]
        p_ = funs_[field_number]
        p_1 = funs_1[field_number]
        field_number += 1

        if p_lagrange:
            q0 = test_functions["NSPFEC"][field_number]
            p0_ = funs_[field_number]
            p0_1 = funs_1[field_number]
            field_number += 1
        else:
            q0 = p0_ = p0_1 = None
    else:
        v = u_ = u_1 = q = q0 = p0_ = p0_1 = None
    if enable_PF:
        psi = test_functions["NSPFEC"][field_number]
        phi_ = funs_[field_number]
        phi_1 = funs_1[field_number]
        field_number += 1
        h = test_functions["NSPFEC"][field_number]
        g_ = funs_[field_number]
        g_1 = funs_1[field_number]
        field_number += 1
    else:
        psi = phi_ = phi_1 = h = g_ = g_1 = 1
    if enable_EC:
        num_solutes = len(test_functions["NSPFEC"])-field_number-1
        b = test_functions["NSPFEC"][field_number:(num_solutes+field_number)]
        c_ = funs_[field_number:(num_solutes+field_number)]
        c_1 = funs_1[field_number:(num_solutes+field_number)]
        U = test_functions["NSPFEC"][num_solutes+field_number]
        V_ = funs_[num_solutes+field_number]
        V_1 = funs_1[num_solutes+field_number]
    else:
        b = c_ = c_1 = U = V_ = V_1 = rho_e_ = 0

    M_ = pf_mobility(phi_, gamma)
    nu_ = ramp(phi_, viscosity)
    veps_ = ramp(phi_, permittivity)
    rho_ = ramp(phi_, density)
    dveps = dramp(permittivity)
    drho = dramp(density)

    dbeta = []  # Diff. in beta
    z = []  # Charge z[species]
    K_ = []  # Diffusivity K[species]
    beta_ = []  # Conc. jump func. beta[species]

    for solute in solutes:
        z.append(solute[1])
        K_.append(ramp(phi_, [solute[2], solute[3]]))
        beta_.append(ramp(phi_, [solute[4], solute[5]]))
        dbeta.append(dramp([solute[4], solute[5]]))

    if enable_EC:
        rho_e_ = sum([c_e*z_e for c_e, z_e in zip(c_, z)])  # Sum of current sol.
        rho_e_1 = sum([c_e*z_e for c_e, z_e in zip(c_1, z)])  # Sum of current sol.

    solver = dict()
    solver["NSPFEC"] = setup_NSPFEC(w_["NSPFEC"], w_1["NSPFEC"],
                                    dirichlet_bcs["NSPFEC"],
                                    neumann_bcs,
                                    boundary_to_mark,
                                    dx, ds, normal,
                                    v, q, q0, psi, h, b, U,
                                    u_, p_, p0_, phi_, g_, c_, V_,
                                    u_1, p_1, p0_1, phi_1, g_1, c_1, V_1,
                                    M_, nu_, veps_, rho_, K_, beta_, rho_e_,
                                    dbeta, dveps, drho,
                                    per_tau, sigma_bar, eps, grav, z,
                                    solutes,
                                    enable_NS, enable_PF, enable_EC,
                                    use_iterative_solvers,
                                    p_lagrange,
                                    q_rhs)
    return dict(solvers=solver)
예제 #2
0
파일: basic.py 프로젝트: mstiegl/BERNAISE
def setup(mesh, test_functions, trial_functions, w_, w_1, ds, dx, normal,
          dirichlet_bcs, neumann_bcs, boundary_to_mark, permittivity, density,
          viscosity, solutes, enable_PF, enable_EC, enable_NS, surface_tension,
          dt, interface_thickness, grav_const, grav_dir, friction_coeff,
          pf_mobility, pf_mobility_coeff, use_iterative_solvers,
          use_pressure_stabilization, comoving_velocity, p_lagrange, q_rhs,
          **namespace):
    """ Set up problem. """
    # Constant
    dim = mesh.geometry().dim()
    sigma_bar = surface_tension * 3. / (2 * math.sqrt(2))
    per_tau = df.Constant(1. / dt)
    grav = df.Constant(tuple(grav_const * np.array(grav_dir[:dim])))
    gamma = pf_mobility_coeff
    eps = interface_thickness
    fric = df.Constant(friction_coeff)
    u_comoving = df.Constant(tuple(comoving_velocity[:dim]))

    # Navier-Stokes
    u_ = p_ = None
    u_1 = p_1 = None
    p0 = q0 = p0_ = p0_1 = None
    if enable_NS:
        u, p = trial_functions["NS"][:2]
        v, q = test_functions["NS"][:2]

        up_ = df.split(w_["NS"])
        up_1 = df.split(w_1["NS"])
        u_, p_ = up_[:2]
        u_1, p_1 = up_1[:2]
        if p_lagrange:
            p0 = trial_functions["NS"][-1]
            q0 = test_functions["NS"][-1]
            p0_ = up_[-1]
            p0_1 = up_1[-1]

    # Phase field
    if enable_PF:
        phi, g = trial_functions["PF"]
        psi, h = test_functions["PF"]

        phi_, g_ = df.split(w_["PF"])
        phi_1, g_1 = df.split(w_1["PF"])
    else:
        # Defaults to phase 1 if phase field is disabled
        phi_ = phi_1 = 1.
        g_ = g_1 = None

    # Electrochemistry
    if enable_EC:
        num_solutes = len(trial_functions["EC"]) - 1
        assert (num_solutes == len(solutes))
        c = trial_functions["EC"][:num_solutes]
        V = trial_functions["EC"][num_solutes]
        b = test_functions["EC"][:num_solutes]
        U = test_functions["EC"][num_solutes]

        cV_ = df.split(w_["EC"])
        cV_1 = df.split(w_1["EC"])
        c_, V_ = cV_[:num_solutes], cV_[num_solutes]
        c_1, V_1 = cV_1[:num_solutes], cV_1[num_solutes]
    else:
        c_ = V_ = c_1 = V_1 = None

    phi_flt_ = unit_interval_filter(phi_)
    phi_flt_1 = unit_interval_filter(phi_1)

    M_ = pf_mobility(phi_flt_, gamma)
    M_1 = pf_mobility(phi_flt_1, gamma)
    mu_ = ramp(phi_flt_, viscosity)
    rho_ = ramp(phi_flt_, density)
    rho_1 = ramp(phi_flt_1, density)
    veps_ = ramp(phi_flt_, permittivity)

    dveps = dramp(permittivity)
    drho = dramp(density)

    dbeta = []  # Diff. in beta
    z = []  # Charge z[species]
    K_ = []  # Diffusivity K[species]
    beta_ = []  # Conc. jump func. beta[species]

    for solute in solutes:
        z.append(solute[1])
        K_.append(ramp(phi_, [solute[2], solute[3]]))
        beta_.append(ramp(phi_, [solute[4], solute[5]]))
        dbeta.append(dramp([solute[4], solute[5]]))

    if enable_EC:
        rho_e = sum([c_e * z_e
                     for c_e, z_e in zip(c, z)])  # Sum of trial func.
        rho_e_ = sum([c_e * z_e
                      for c_e, z_e in zip(c_, z)])  # Sum of curr. sol.
    else:
        rho_e_ = None

    solvers = dict()
    if enable_PF:
        solvers["PF"] = setup_PF(w_["PF"], phi, g, psi, h, dx, ds, normal,
                                 dirichlet_bcs["PF"], neumann_bcs,
                                 boundary_to_mark, phi_1, u_1, M_1, c_1, V_1,
                                 per_tau, sigma_bar, eps, dbeta, dveps,
                                 enable_NS, enable_EC, use_iterative_solvers,
                                 q_rhs)

    if enable_EC:
        solvers["EC"] = setup_EC(w_["EC"], c, V, b, U, rho_e, dx, ds, normal,
                                 dirichlet_bcs["EC"], neumann_bcs,
                                 boundary_to_mark, c_1, u_1, K_, veps_,
                                 phi_flt_, solutes, per_tau, z, dbeta,
                                 enable_NS, enable_PF, use_iterative_solvers,
                                 q_rhs)

    if enable_NS:
        solvers["NS"] = setup_NS(w_["NS"], u, p, v, q, p0, q0, dx, ds, normal,
                                 dirichlet_bcs["NS"], neumann_bcs,
                                 boundary_to_mark, u_1, phi_flt_, rho_, rho_1,
                                 g_, M_, mu_, rho_e_, c_, V_, c_1, V_1, dbeta,
                                 solutes, per_tau, drho, sigma_bar, eps, dveps,
                                 grav, fric, u_comoving, enable_PF, enable_EC,
                                 use_iterative_solvers,
                                 use_pressure_stabilization, p_lagrange, q_rhs)
    return dict(solvers=solvers)
예제 #3
0
파일: TDLUES.py 프로젝트: mstiegl/BERNAISE
def unpack_quantities(surface_tension, grav_const, grav_dir, pf_mobility_coeff,
                      pf_mobility, interface_thickness, density, viscosity,
                      permittivity, trial_functions, test_functions, w_, w_1,
                      solutes, dt, enable_EC, enable_PF, enable_NS):
    """ """
    # Constant
    sigma_bar = surface_tension * 3. / (2 * math.sqrt(2))
    per_tau = df.Constant(1. / dt)
    grav = df.Constant(tuple(grav_const * np.array(grav_dir)))
    gamma = pf_mobility_coeff
    eps = interface_thickness
    rho_min = min(density) if enable_PF else density[0]

    # Navier-Stokes
    if enable_NS:
        u = trial_functions["NSu"]
        p = trial_functions["NSp"]
        v = test_functions["NSu"]
        q = test_functions["NSp"]

        u_ = w_["NSu"]
        p_ = w_["NSp"]
        u_1 = w_1["NSu"]
        p_1 = w_1["NSp"]
    else:
        u = p = v = q = None
        u_ = p_ = None
        u_1 = p_1 = None

    # Phase field
    if enable_PF:
        phi, g = trial_functions["PF"]
        psi, h = test_functions["PF"]

        phi_, g_ = df.split(w_["PF"])
        phi_1, g_1 = df.split(w_1["PF"])
    else:
        # Defaults to phase 1 if phase field is disabled
        phi = 1.
        g = psi = h = None
        phi_ = phi_1 = 1.
        g_ = g_1 = None

    # Electrochemistry
    if enable_EC:
        num_solutes = len(trial_functions["EC"]) - 1
        assert (num_solutes == len(solutes))
        c = trial_functions["EC"][:num_solutes]
        V = trial_functions["EC"][num_solutes]
        b = test_functions["EC"][:num_solutes]
        U = test_functions["EC"][num_solutes]

        cV_ = df.split(w_["EC"])
        cV_1 = df.split(w_1["EC"])
        c_, V_ = cV_[:num_solutes], cV_[num_solutes]
        c_1, V_1 = cV_1[:num_solutes], cV_1[num_solutes]
    else:
        c = V = b = U = c_ = V_ = c_1 = V_1 = None

    phi_flt_ = unit_interval_filter(phi_)
    phi_flt_1 = unit_interval_filter(phi_1)
    # phi_flt_ = phi_
    # phi_flt_1 = phi_1
    #phi_ = phi_flt_
    phi_1 = phi_flt_1

    M_ = pf_mobility(phi_flt_, gamma)
    M_1 = pf_mobility(phi_flt_1, gamma)
    nu_ = ramp(phi_flt_, viscosity)
    rho_ = ramp(phi_flt_, density)
    veps_ = ramp(phi_flt_, permittivity)

    rho_1 = ramp(phi_flt_1, density)
    nu_1 = ramp(phi_flt_1, viscosity)

    dveps = dramp(permittivity)
    drho = dramp(density)

    dbeta = []  # Diff. in beta
    z = []  # Charge z[species]
    K_ = []  # Diffusivity K[species]
    beta_ = []  # Conc. jump func. beta[species]

    for solute in solutes:
        z.append(solute[1])
        K_.append(ramp_geometric(phi_flt_, [solute[2], solute[3]]))
        beta_.append(ramp(phi_flt_, [solute[4], solute[5]]))
        dbeta.append(dramp([solute[4], solute[5]]))

    if enable_EC:
        rho_e = sum([c_e * z_e for c_e, z_e in zip(c, z)])
        rho_e_ = sum([c_e_ * z_e for c_e_, z_e in zip(c_, z)])
        rho_e_1 = sum([c_e_1 * z_e for c_e_1, z_e in zip(c_1, z)])
    else:
        rho_e = rho_e_ = rho_e_1 = None

    return (sigma_bar, per_tau, grav, gamma, eps, rho_min, u, p, v, q, u_, p_,
            u_1, p_1, phi, g, psi, h, phi_, g_, phi_1, g_1, c, V, b, U, c_, V_,
            c_1, V_1, phi_flt_, phi_flt_1, M_, M_1, nu_, rho_, veps_, rho_1,
            nu_1, dveps, drho, dbeta, z, K_, beta_, rho_e, rho_e_, rho_e_1)
예제 #4
0
def setup(tstep, test_functions, trial_functions, w_, w_1, ds, dx, normal,
          dirichlet_bcs, neumann_bcs, boundary_to_mark, permittivity, density,
          viscosity, solutes, enable_PF, enable_EC, enable_NS, surface_tension,
          dt, interface_thickness, grav_const, pf_mobility, pf_mobility_coeff,
          use_iterative_solvers, use_pressure_stabilization, q_rhs,
          **namespace):
    """ Set up problem. """
    # Constant
    sigma_bar = surface_tension * 3. / (2 * math.sqrt(2))
    per_tau = df.Constant(1. / dt)
    grav = df.Constant((0., -grav_const))
    gamma = pf_mobility_coeff
    eps = interface_thickness

    # Navier-Stokes
    if enable_NS:
        u = trial_functions["NSu"]
        p = trial_functions["NSp"]
        v = test_functions["NSu"]
        q = test_functions["NSp"]

        u_ = w_["NSu"]
        p_ = w_["NSp"]
        u_1 = w_1["NSu"]
        p_1 = w_1["NSp"]

    # Phase field
    if enable_PF:
        phi, g = trial_functions["PF"]
        psi, h = test_functions["PF"]

        phi_, g_ = df.split(w_["PF"])
        phi_1, g_1 = df.split(w_1["PF"])
    else:
        # Defaults to phase 1 if phase field is disabled
        phi_ = phi_1 = 1.

    # Electrochemistry
    if enable_EC:
        num_solutes = len(trial_functions["EC"]) - 1
        assert (num_solutes == len(solutes))
        c = trial_functions["EC"][:num_solutes]
        V = trial_functions["EC"][num_solutes]
        b = test_functions["EC"][:num_solutes]
        U = test_functions["EC"][num_solutes]

        cV_ = df.split(w_["EC"])
        cV_1 = df.split(w_1["EC"])
        c_, V_ = cV_[:num_solutes], cV_[num_solutes]
        c_1, V_1 = cV_1[:num_solutes], cV_1[num_solutes]
    else:
        c_, V_ = None, None
        c_1, V_1 = None, None

    M_ = pf_mobility(unit_interval_filter(phi_), gamma)
    M_1 = pf_mobility(unit_interval_filter(phi_1), gamma)
    nu_ = ramp(unit_interval_filter(phi_), viscosity)
    rho_ = ramp(unit_interval_filter(phi_), density)
    veps_ = ramp(unit_interval_filter(phi_), permittivity)

    rho_1 = ramp(unit_interval_filter(phi_1), density)

    dveps = dramp(permittivity)
    drho = dramp(density)

    dbeta = []  # Diff. in beta
    z = []  # Charge z[species]
    K_ = []  # Diffusivity K[species]
    beta_ = []  # Conc. jump func. beta[species]

    for solute in solutes:
        z.append(solute[1])
        K_.append(ramp(phi_, [solute[2], solute[3]]))
        beta_.append(ramp(phi_, [solute[4], solute[5]]))
        dbeta.append(dramp([solute[4], solute[5]]))

    if enable_EC:
        rho_e = sum([c_e * z_e
                     for c_e, z_e in zip(c, z)])  # Sum of trial functions
        rho_e_ = sum([c_e * z_e
                      for c_e, z_e in zip(c_, z)])  # Sum of current sol.
    else:
        rho_e = None
        rho_e_ = None

    if tstep == 0 and enable_NS:
        solve_initial_pressure(w_["NSp"], p, q, u, v, dirichlet_bcs["NSp"], M_,
                               g_, phi_, rho_, rho_e_, V_, drho, sigma_bar,
                               eps, grav, dveps, enable_PF, enable_EC)

    solvers = dict()
    if enable_PF:
        # solvers["PF"] = setup_PF(w_["PF"], phi, g, psi, h, dirichlet_bcs["PF"],
        #                          phi_1, u_1, M_1, c_1, V_1,
        #                          per_tau, sigma_bar, eps, dbeta, dveps,
        #                          enable_NS, enable_EC,
        #                          use_iterative_solvers)
        solvers["PF"] = setup_PF(w_["PF"], phi, g, psi, h, dx, ds, normal,
                                 dirichlet_bcs["PF"], neumann_bcs,
                                 boundary_to_mark, phi_1, u_1, M_1, c_1, V_1,
                                 per_tau, sigma_bar, eps, dbeta, dveps,
                                 enable_NS, enable_EC, use_iterative_solvers,
                                 q_rhs)

    if enable_EC:
        solvers["EC"] = setup_EC(w_["EC"], c, V, b, U, rho_e,
                                 dirichlet_bcs["EC"], c_1, u_1, K_, veps_,
                                 phi_, per_tau, z, dbeta, enable_NS, enable_PF,
                                 use_iterative_solvers)

    if enable_NS:
        solvers["NSu"] = setup_NSu(u, v, u_, p_, dirichlet_bcs["NSu"], u_1,
                                   p_1, phi_, rho_, rho_1, g_, M_, nu_, rho_e_,
                                   V_, dt, drho, sigma_bar, eps, dveps, grav,
                                   enable_PF, enable_EC)
        solvers["NSp"] = setup_NSp(p, q, dirichlet_bcs["NSp"], u_, p_, p_1,
                                   rho_, dt)

    return dict(solvers=solvers)
예제 #5
0
def setup(test_functions, trial_functions, w_, w_1, bcs, permittivity, density,
          viscosity, solutes, enable_PF, enable_EC, enable_NS, surface_tension,
          dt, interface_thickness, grav_const, pf_mobility_coeff, pf_mobility,
          use_iterative_solvers, **namespace):
    """ Set up problem. """
    # Constant
    sigma_bar = surface_tension * 3. / (2 * math.sqrt(2))
    per_tau = df.Constant(1. / dt)
    grav = df.Constant((0., -grav_const))
    gamma = pf_mobility_coeff
    eps = interface_thickness

    funs_ = df.split(w_["NSPFEC"])
    funs_1 = df.split(w_1["NSPFEC"])
    field_number = 0
    if enable_NS:
        v = test_functions["NSPFEC"][field_number]
        u_ = funs_[field_number]
        u_1 = funs_1[field_number]
        field_number += 1
        q = test_functions["NSPFEC"][field_number]
        p_ = funs_[field_number]
        p_1 = funs_1[field_number]
        field_number += 1
    if enable_PF:
        psi = test_functions["NSPFEC"][field_number]
        phi_ = funs_[field_number]
        phi_1 = funs_1[field_number]
        field_number += 1
        h = test_functions["NSPFEC"][field_number]
        g_ = funs_[field_number]
        g_1 = funs_1[field_number]
        field_number += 1
    if enable_EC:
        num_solutes = len(test_functions["NSPFEC"]) - field_number - 1
        b = test_functions["NSPFEC"][field_number:(num_solutes + field_number)]
        c_ = funs_[field_number:(num_solutes + field_number)]
        c_1 = funs_1[field_number:(num_solutes + field_number)]
        U = test_functions["NSPFEC"][num_solutes + field_number]
        V_ = funs_[num_solutes + field_number]
        V_1 = funs_1[num_solutes + field_number]

    M_ = pf_mobility(phi_, gamma)
    M_1 = pf_mobility(phi_1, gamma)
    nu_ = ramp(phi_, viscosity)
    nu_1 = ramp(phi_1, viscosity)
    veps_ = ramp(phi_, permittivity)
    veps_1 = ramp(phi_1, permittivity)
    rho_ = ramp(phi_, density)
    rho_1 = ramp(phi_1, density)
    dveps = dramp(permittivity)
    drho = dramp(density)

    dbeta = []  # Diff. in beta
    z = []  # Charge z[species]
    K_ = []  # Diffusivity K[species]
    K_1 = []
    beta_ = []  # Conc. jump func. beta[species]
    beta_1 = []

    for solute in solutes:
        z.append(solute[1])
        K_.append(ramp(phi_, [solute[2], solute[3]]))
        K_1.append(ramp(phi_1, [solute[2], solute[3]]))
        beta_.append(ramp(phi_, [solute[4], solute[5]]))
        beta_1.append(ramp(phi_, [solute[4], solute[5]]))
        dbeta.append(dramp([solute[4], solute[5]]))

    if enable_EC:
        rho_e_ = sum([c_e * z_e
                      for c_e, z_e in zip(c_, z)])  # Sum of curr. sol.
        rho_e_1 = sum([c_e * z_e for c_e, z_e in zip(c_1, z)])  # prev sol.

    solver = dict()
    solver["NSPFEC"] = setup_NSPFEC(
        w_["NSPFEC"], w_1["NSPFEC"], bcs["NSPFEC"], trial_functions["NSPFEC"],
        v, q, psi, h, b, U, u_, p_, phi_, g_, c_, V_, u_1, p_1, phi_1, g_1,
        c_1, V_1, M_, nu_, veps_, rho_, K_, rho_e_, M_1, nu_1, veps_1, rho_1,
        K_1, rho_e_1, dbeta, dveps, drho, per_tau, sigma_bar, eps, grav, z,
        enable_NS, enable_PF, enable_EC, use_iterative_solvers)
    return dict(solvers=solver)
예제 #6
0
def method(ts, dt=0, extra_boundaries="", **kwargs):
    """ Plot flux in time. """
    info_cyan("Plot flux in time.")

    params = ts.get_parameters()
    steps = get_steps(ts, dt)

    problem = params["problem"]
    info("Problem: {}".format(problem))

    boundary_to_mark, ds = fetch_boundaries(ts, problem, params,
                                            extra_boundaries)

    x_ = ts.functions()

    if params["enable_NS"]:
        u = x_["u"]
    else:
        u = df.Constant(0.)

    if params["enable_PF"]:
        phi = x_["phi"]
        g = x_["g"]
        exec("from problems.{} import pf_mobility".format(problem))
        M = pf_mobility(phi, params["pf_mobility_coeff"])
    else:
        phi = 1.
        g = df.Constant(0.)
        M = df.Constant(0.)

    solutes = params["solutes"]
    c = []
    c_grad_g_c = []
    if params["enable_EC"]:
        V = x_["V"]
    else:
        V = df.Constant(0.)

    dbeta = []  # Diff. in beta
    z = []  # Charge z[species]
    K = []  # Diffusivity K[species]
    beta = []  # Conc. jump func. beta[species]

    for solute in solutes:
        ci = x_[solute[0]]
        dbetai = dramp([solute[4], solute[5]])
        c.append(ci)
        z.append(solute[1])
        K.append(ramp(phi, [solute[2], solute[3]]))
        beta.append(ramp(phi, [solute[4], solute[5]]))
        dbeta.append(dbetai)
        # THIS HAS NOT BEEN GENERALIZED!
        c_grad_g_ci = df.grad(ci) + solute[1] * ci * df.grad(V)
        if params["enable_PF"]:
            c_grad_g_ci += dbetai * df.grad(phi)
        c_grad_g_c.append(c_grad_g_ci)

    nu = ramp(phi, params["viscosity"])
    veps = ramp(phi, params["permittivity"])
    rho = ramp(phi, params["density"])

    dveps = dramp(params["permittivity"])
    drho = dramp(params["density"])

    t = np.zeros(len(steps))

    # Define the fluxes
    fluxes = dict()
    fluxes["Velocity"] = u
    fluxes["Phase"] = phi * u
    fluxes["Mass"] = rho * x_["u"]
    if params["enable_PF"]:
        fluxes["Phase"] += -M * df.grad(g)
        fluxes["Mass"] += -drho * M * df.grad(g)

    if params["enable_EC"]:
        for i, solute in enumerate(solutes):
            fluxes["Solute {}".format(solute[0])] = K[i] * c_grad_g_c[i]
        fluxes["E-field"] = -df.grad(V)

    data = dict()
    for boundary_name in boundary_to_mark:
        data[boundary_name] = dict()
        for flux_name in fluxes:
            data[boundary_name][flux_name] = np.zeros(len(steps))

    n = df.FacetNormal(ts.mesh)

    for i, step in enumerate(steps):
        info("Step {} of {}".format(step, len(ts)))

        for field in x_:
            ts.update(x_[field], field, step)

        for boundary_name, (mark, k) in boundary_to_mark.items():
            for flux_name, flux in fluxes.items():
                data[boundary_name][flux_name][i] = df.assemble(
                    df.dot(flux, n) * ds[k](mark))

        t[i] = ts.times[step]

    savedata = dict()
    flux_keys = sorted(fluxes.keys())
    for boundary_name in boundary_to_mark:
        savedata[boundary_name] = np.array(
            list(
                zip(
                    steps, t, *[
                        data[boundary_name][flux_name]
                        for flux_name in flux_keys
                    ])))

    if rank == 0:
        header = "Step\tTime\t" + "\t".join(flux_keys)
        for boundary_name in boundary_to_mark:
            with open(
                    os.path.join(ts.analysis_folder,
                                 "flux_in_time_{}.dat".format(boundary_name)),
                    "w") as outfile:
                np.savetxt(outfile, savedata[boundary_name], header=header)