Exemple #1
0
class nozzle(euler1d):
    """
    Class nozzle for euler equations with section term -1/A dA/dx (rho u, rho u2, rho u Ht)

    attributes:

    """
    _bcdict = base.methoddict('bc_')
    _vardict = base.methoddict()
    _numfluxdict = base.methoddict('numflux_')

    def __init__(self, sectionlaw, gamma=1.4, source=None):
        nozsrc = [self.src_mass, self.src_mom, self.src_energy]
        allsrc = nozsrc  # init all sources to nozzle sources
        if source:  # additional sources ?
            for i, isrc in enumerate(source):
                if isrc:
                    allsrc[i] = lambda x, q: isrc(x, q) + nozsrc[i](x, q)
        euler1d.__init__(self, gamma=gamma, source=allsrc)
        self.sectionlaw = sectionlaw
        self._bcdict.merge(nozzle._bcdict)
        self._vardict.merge(nozzle._vardict)
        self._numfluxdict.merge(nozzle._numfluxdict)

    def initdisc(self, mesh):
        self._xc = mesh.centers()
        self.geomterm = 1./self.sectionlaw(self._xc)* \
                    (self.sectionlaw(mesh.xf[1:mesh.ncell+1])-self.sectionlaw(mesh.xf[0:mesh.ncell])) / \
                    (mesh.xf[1:mesh.ncell+1]-mesh.xf[0:mesh.ncell])
        return

    @_vardict.register()
    def massflow(self, qdata):
        return qdata[1] * self.sectionlaw(self._xc)

    def src_mass(self, x, qdata):
        return -self.geomterm * qdata[1]

    def src_mom(self, x, qdata):
        return -self.geomterm * qdata[1]**2 / qdata[0]

    def src_energy(self, x, qdata):
        ec = 0.5 * qdata[1]**2 / qdata[0]
        return -self.geomterm * qdata[1] * (
            (qdata[2] - ec) * self.gamma + ec) / qdata[0]
Exemple #2
0
class a():
    _dico = methoddict('f_')

    def __init__(self):
        self._dico = a._dico.copy()

    @_dico.register()
    def f_ma(self):
        return

    @_dico.register('')
    def m(self):
        return
Exemple #3
0
class b(a):
    _dico = methoddict()

    def __init__(self):
        a.__init__(self)
        self._dico.merge(b._dico)

    @_dico.register('f_')
    def f_mb(self):
        return

    @_dico.register()
    def m(self):
        return
Exemple #4
0
class model(base.model):
    """
    Class model for convection

    attributes:

    """
    _vardict = base.methoddict()

    def __init__(self, convcoef):
        base.model.__init__(self, name='convection', neq=1)
        self.has_firstorder_terms = 1
        self.convcoef = convcoef
        self.islinear = 1
        self.shape = [1]
        self._vardict.merge(model._vardict)

    def cons2prim(self, qdata):
        return [1 * d for d in qdata]

    def prim2cons(self, pdata):
        return [1 * d for d in pdata]

    def numflux(self, name, pL, pR):
        """
        >>> model(2.).numflux([10.],[50.])
        [20.0]
        >>> model(-2.).numflux([10.],[50.])
        [-100.0]
        """
        return [
            self.convcoef * (pL[0] + pR[0]) / 2. - abs(self.convcoef) *
            (pR[0] - pL[0]) / 2.
        ]

    def timestep(self, pdata, dx, condition):
        "computation of timestep: data is not used, dx is an array of cell sizes, condition is the CFL number"
        return condition * dx / abs(self.convcoef)

    @_vardict.register()
    def q(self, qdata):
        return qdata[0].copy()
Exemple #5
0
class c(a):
    _dico = methoddict('f_')

    def __init__(self):
        a.__init__(self)
        self._dico.merge(c._dico)

    @_dico.register(name='zmc')
    def f_mc(self):
        return

    try:

        @_dico.register()
        def m(self):
            return
    except LookupError as error:
        logging.traceback.print_exc()
        print('Error is ignored')

        @_dico.register(name='zm')
        def m(self):
            return
Exemple #6
0
class euler2d(euler):
    """
    Class model for 2D euler equations
    """
    _bcdict = base.methoddict('bc_')
    _vardict = base.methoddict()
    _numfluxdict = base.methoddict('numflux_')

    def __init__(self, gamma=1.4, source=None):
        euler.__init__(self, gamma=gamma, source=source)
        self.shape = [1, 2, 1]
        self._bcdict.merge(euler2d._bcdict)
        self._vardict.merge(euler2d._vardict)
        self._numfluxdict.merge(euler2d._numfluxdict)

    def _derived_fromprim(self, pdata, dir):
        """
        returns rho, un, V, p, H, c2
        """
        c2 = self.gamma * pdata[2] / pdata[0]
        un = _vec_dot_vec(pdata[1], dir)
        H = c2 / (self.gamma - 1.) + .5 * _vecsqrmag(pdata[1])
        return pdata[0], un, pdata[1], pdata[2], H, c2

    @_vardict.register()
    def velocity_x(self, qdata):  # returns (rho ux)/rho
        return qdata[1][0, :] / qdata[0]

    @_vardict.register()
    def velocity_y(self, qdata):  # returns (rho uy)/rho
        return qdata[1][1, :] / qdata[0]

    @_vardict.register()
    def mach(self, qdata):
        rhoUmag = _vecmag(qdata[1])
        return rhoUmag / np.sqrt(self.gamma *
                                 ((self.gamma - 1.0) *
                                  (qdata[0] * qdata[2] - 0.5 * rhoUmag**2)))

    def _Roe_average(self, rhoL, unL, UL, HL, rhoR, unR, UR, HR):
        """returns Roe averaged rho, u, usound"""
        # Roe's averaging
        Rrho = np.sqrt(rhoR / rhoL)
        tmp = 1.0 / (1.0 + Rrho)
        unRoe = tmp * (unL + unR * Rrho)
        URoe = tmp * (UL + UR * Rrho)
        hRoe = tmp * (HL + HR * Rrho)
        cRoe = np.sqrt((hRoe - 0.5 * _vecsqrmag(URoe)) * (self.gamma - 1.))
        return Rrho, unRoe, cRoe

    @_numfluxdict.register(name='centered')
    @_numfluxdict.register()
    def numflux_centeredflux(self, pdataL, pdataR,
                             dir):  # centered flux ; pL[ieq][face]
        rhoL, unL, VL, pL, HL, cL2 = self._derived_fromprim(pdataL, dir)
        rhoR, unR, VR, pR, HR, cR2 = self._derived_fromprim(pdataR, dir)
        # final flux
        Frho = .5 * (rhoL * unL + rhoR * unR)
        Frhou = .5 * ((rhoL * unL) * VL + pL * dir +
                      (rhoR * unR) * VR + pR * dir)
        FrhoE = .5 * ((rhoL * unL * HL) + (rhoR * unR * HR))
        return [Frho, Frhou, FrhoE]

    @_numfluxdict.register(name='hlle')
    def numflux_hlle(self, pdataL, pdataR,
                     dir):  # HLLE Riemann solver ; pL[ieq][face]
        rhoL, unL, VL, pL, HL, cL2 = self._derived_fromprim(pdataL, dir)
        rhoR, unR, VR, pR, HR, cR2 = self._derived_fromprim(pdataR, dir)
        # The HLLE Riemann solver
        etL = HL - pL / rhoL
        etR = HR - pR / rhoR
        # Roe's averaging
        Rrho, uRoe, cRoe = self._Roe_average(rhoL, unL, VL, HL, rhoR, unR, VR,
                                             HR)
        # max HLL 2 waves "velocities"
        sL = np.minimum(0., np.minimum(uRoe - cRoe, unL - np.sqrt(cL2)))
        sR = np.maximum(0., np.maximum(uRoe + cRoe, unR + np.sqrt(cR2)))
        # final flux
        Frho = (sR * rhoL * unL - sL * rhoR * unR + sL * sR *
                (rhoR - rhoL)) / (sR - sL)
        Frhou = (sR * ((rhoL * unL) * VL + pL * dir) - sL *
                 ((rhoR * unR) * VR + pR * dir) + sL * sR *
                 (rhoR * VR - rhoL * VL)) / (sR - sL)
        FrhoE = (sR * (rhoL * unL * HL) - sL * (rhoR * unR * HR) + sL * sR *
                 (rhoR * etR - rhoL * etL)) / (sR - sL)
        return [Frho, Frhou, FrhoE]

    @_bcdict.register()
    def bc_sym(self, dir, data, param):
        "symmetry boundary condition, for inviscid equations, it is equivalent to a wall, do not need user parameters"
        VL = data[1]
        Vn = _vec_dot_vec(VL, dir)
        VR = VL - 2.0 * (Vn * dir)
        return [data[0], VR, data[2]]

    @_bcdict.register()
    def bc_insub(self, dir, data, param):
        #needed parameters : ptot, rttot
        g = self.gamma
        gmu = g - 1.
        p = data[2]
        m2 = np.maximum(0., ((param['ptot'] / p)**(gmu / g) - 1.) * 2. / gmu)
        rh = param['ptot'] / param['rttot'] / (1. + .5 * gmu * m2)**(1. / gmu)
        return [rh, _sca_mult_vec(-np.sqrt(g * p * m2 / rh), dir), p]

    @_bcdict.register()
    def bc_insup(self, dir, data, param):
        # needed parameters : ptot, rttot
        g = self.gamma
        gmu = g - 1.
        p = param['p']
        if 'angle' in param:
            ang = np.deg2rad(param['angle'])
            dir_in = np.full_like(dir, [[np.cos(ang)], [np.sin(ang)]])
        else:
            dir_in = -dir
        m2 = np.maximum(0., ((param['ptot'] / p)**(gmu / g) - 1.) * 2. / gmu)
        rh = param['ptot'] / param['rttot'] / (1. + .5 * gmu * m2)**(1. / gmu)
        return [rh, _sca_mult_vec(np.sqrt(g * p * m2 / rh), dir_in), p]

    @_bcdict.register()
    def bc_outsub(self, dir, data, param):
        return [data[0], data[1], param['p']]

    @_bcdict.register()
    def bc_outsup(self, dir, data, param):
        return data
Exemple #7
0
class euler1d(euler):
    """
    Class model for 1D euler equations
    """
    _bcdict = base.methoddict('bc_')
    _vardict = base.methoddict()
    _numfluxdict = base.methoddict('numflux_')

    def __init__(self, gamma=1.4, source=None):
        euler.__init__(self, gamma=gamma, source=source)
        self.shape = [1, 1, 1]
        self._bcdict.merge(euler1d._bcdict)
        self._vardict.merge(euler1d._vardict)
        self._numfluxdict.merge(euler1d._numfluxdict)

    def _derived_fromprim(self, pdata, dir):
        """
        returns rho, un, V, c2, H
        'dir' is ignored
        """
        c2 = self.gamma * pdata[2] / pdata[0]
        H = c2 / (self.gamma - 1.) + .5 * pdata[1]**2
        return pdata[0], pdata[1], pdata[1], c2, H

    @_vardict.register()
    def massflow(self, qdata):  # for 1D model only
        return qdata[1].copy()

    @_bcdict.register()
    def bc_sym(self, dir, data, param):
        "symmetry boundary condition, for inviscid equations, it is equivalent to a wall, do not need user parameters"
        return [data[0], -data[1], data[2]]

    @_bcdict.register()
    def bc_insub(self, dir, data, param):
        g = self.gamma
        gmu = g - 1.
        p = data[2]
        m2 = np.maximum(0., ((param['ptot'] / p)**(gmu / g) - 1.) * 2. / gmu)
        rh = param['ptot'] / param['rttot'] / (1. + .5 * gmu * m2)**(1. / gmu)
        return [rh, -dir * np.sqrt(g * m2 * p / rh), p]

    @_bcdict.register()
    def bc_insub_cbc(self, dir, data, param):
        g = self.gamma
        gmu = g - 1.
        p = data[2]
        invcm = data[1] + dir * 2 * np.sqrt(g * p / data[0]) / gmu
        adiscri = g * (g + 1) / gmu * param['rttot'] - .5 * gmu * invcm**2
        a1 = (dir * invcm + np.sqrt(adiscri)) * gmu / (g + 1.)
        u1 = invcm - dir * 2 * a1 / gmu
        f_m1sqr = 1. + .5 * gmu * (u1 / a1)**2
        rh1 = param['ptot'] / param['rttot'] / f_m1sqr**(1. / gmu)
        p1 = param['ptot'] / f_m1sqr**(g / gmu)
        return [rh1, u1, p1]

    @_bcdict.register()
    def bc_insup(self, dir, data, param):
        # expected parameters are 'ptot', 'rttot' and 'p'
        g = self.gamma
        gmu = g - 1.
        p = param['p']
        m2 = np.maximum(0., ((param['ptot'] / p)**(gmu / g) - 1.) * 2. / gmu)
        rh = param['ptot'] / param['rttot'] / (1. + .5 * gmu * m2)**(1. / gmu)
        return [rh, -dir * np.sqrt(g * m2 * p / rh), p]

    @_bcdict.register(name='outsub')
    @_bcdict.register()
    def bc_outsub_prim(self, dir, data, param):
        return [data[0], data[1], param['p']]

    @_bcdict.register()
    def bc_outsub_qtot(self, dir, data, param):
        g = self.gamma
        gmu = g - 1.
        m2 = data[1]**2 / (g * data[2] / data[0])
        fm2 = 1. + .5 * gmu * m2
        rttot = data[2] / data[0] * fm2
        ptot = data[2] * fm2**(g / gmu)
        # right (external) state
        p = param['p']
        m2 = np.maximum(0., ((ptot / p)**(gmu / g) - 1.) * 2. / gmu)
        rho = ptot / rttot / (1. + .5 * gmu * m2)**(1. / gmu)
        return [rho, dir * np.sqrt(g * m2 * p / rho), p]

    @_bcdict.register()
    def bc_outsub_rh(self, dir, data, param):
        g = self.gamma
        gmu = g - 1.
        # pratio > Ms > Ws/a0
        p0 = data[2]
        pratio = param['p'] / p0
        u0 = data[1]
        # relative shock Mach number Ms=(u0-Ws)/a0
        Ms2 = 1. + (pratio - 1.) * (g + 1.) / (2. * g)
        rhoratio = ((g + 1.) * Ms2) / (2. + gmu * Ms2)
        Ws = u0 - dir * np.sqrt(g * p0 / data[0] * Ms2)
        # right (external) state
        p1 = param['p']
        u1 = Ws + (u0 - Ws) / rhoratio
        rho1 = data[0] * rhoratio
        return [rho1, u1, p1]

    @_bcdict.register()
    def bc_outsub_nrcbc(self, dir, data, param):
        g = self.gamma
        gmu = g - 1.
        # 0 and 1 stand for internal/external
        p1 = param['p']
        # isentropic invariant p/rho**gam = cst
        rho1 = data[0] * (p1 / data[2])**(1. / g)
        # C- invariant (or C+ according to dir)
        a0 = np.sqrt(g * data[2] / data[0])
        a1 = np.sqrt(g * p1 / rho1)
        u1 = data[1] + dir * 2 / gmu * (a1 - a0)
        return [rho1, u1, p1]

    @_bcdict.register()
    def bc_outsup(self, dir, data, param):
        return data
Exemple #8
0
class euler(base.model):
    """
    Class model for euler equations

    attributes:

    """
    _bcdict = base.methoddict(
        'bc_')  # dict and associated decorator method to register BC
    _vardict = base.methoddict()
    _numfluxdict = base.methoddict('numflux_')

    def __init__(self, gamma=1.4, source=None):
        base.model.__init__(self, name='euler', neq=3)
        self.islinear = 0
        self.shape = [1, 1, 1]
        self.gamma = gamma
        self.source = source
        self._bcdict.merge(euler._bcdict)
        self._vardict.merge(euler._vardict)
        self._numfluxdict.merge(euler._numfluxdict)

    def cons2prim(self, qdata):  # qdata[ieq][cell] :
        """
        Primitives variables are rho, u, p
        >>> model().cons2prim([[5.], [10.], [20.]])
        True
        """
        rho = qdata[0]
        u = qdata[1] / qdata[0]
        p = self.pressure(qdata)
        pdata = [rho, u, p]
        return pdata

    def prim2cons(self, pdata):  # qdata[ieq][cell] :
        """
        >>> model().prim2cons([[2.], [4.], [10.]]) == [[2.], [8.], [41.]]
        True
        """
        V2 = pdata[1]**2 if pdata[1].ndim == 1 else _vecsqrmag(pdata[1])
        rhoe = pdata[2] / (self.gamma - 1.) + .5 * pdata[0] * V2
        return [pdata[0], _sca_mult_vec(pdata[0], pdata[1]), rhoe]

    @_vardict.register()
    def density(self, qdata):
        return qdata[0].copy()

    @_vardict.register()
    def pressure(self,
                 qdata):  # returns (gam-1)*( rho.et) - .5 * (rho.u)**2 / rho )
        return (self.gamma - 1.0) * (qdata[2] - self.kinetic_energy(qdata))

    @_vardict.register()
    def velocity(
            self,
            qdata):  # returns (rho u)/rho, works for both scalar and vector
        return qdata[1] / qdata[0]

    @_vardict.register()
    def velocitymag(
            self,
            qdata):  # returns mag(rho u)/rho, depending if scalar or vector
        return np.abs(qdata[1]) / qdata[0] if qdata[1].ndim == 1 else _vecmag(
            qdata[1]) / qdata[0]

    @_vardict.register(name="kinetic-energy")
    @_vardict.register()
    def kinetic_energy(self, qdata):
        """volumic kinetic energy"""
        return .5 * qdata[1]**2 / qdata[0] if qdata[
            1].ndim == 1 else .5 * _vecsqrmag(qdata[1]) / qdata[0]

    @_vardict.register()
    def asound(self, qdata):
        return np.sqrt(self.gamma * self.pressure(qdata) / qdata[0])

    @_vardict.register()
    def mach(self, qdata):
        return qdata[1] / np.sqrt(self.gamma *
                                  ((self.gamma - 1.0) *
                                   (qdata[0] * qdata[2] - 0.5 * qdata[1]**2)))

    @_vardict.register()
    def entropy(self, qdata):  # S/r
        return np.log(
            self.pressure(qdata) / qdata[0]**self.gamma) / (self.gamma - 1.)

    @_vardict.register()
    def enthalpy(self, qdata):
        return (qdata[2] -
                0.5 * qdata[1]**2 / qdata[0]) * self.gamma / qdata[0]

    @_vardict.register()
    def ptot(self, qdata):
        gm1 = self.gamma - 1.
        return self.pressure(qdata) * (1. + .5 * gm1 * self.mach(qdata)**2)**(
            self.gamma / gm1)

    @_vardict.register()
    def rttot(self, qdata):
        ec = self.kinetic_energy(qdata)
        return ((qdata[2] - ec) * self.gamma +
                ec) / qdata[0] / self.gamma * (self.gamma - 1.)

    @_vardict.register()
    def htot(self, qdata):
        ec = self.kinetic_energy(qdata)
        return ((qdata[2] - ec) * self.gamma + ec) / qdata[0]

    def _Roe_average(self, rhoL, uL, HL, rhoR, uR, HR):
        """returns Roe averaged rho, u, usound"""
        # Roe's averaging
        Rrho = np.sqrt(rhoR / rhoL)
        tmp = 1.0 / (1.0 + Rrho)
        #velRoe = tmp*(uL + uR*Rrho)
        uRoe = tmp * (uL + uR * Rrho)
        hRoe = tmp * (HL + HR * Rrho)
        cRoe = np.sqrt((hRoe - 0.5 * uRoe**2) * (self.gamma - 1.))
        return Rrho, uRoe, cRoe

    def numflux(self, name, pdataL, pdataR, dir=None):
        if name is None: name = 'hllc'
        return (self._numfluxdict.dict[name])(self, pdataL, pdataR, dir)

    @_numfluxdict.register(name='centered')
    @_numfluxdict.register()
    def numflux_centeredflux(self,
                             pdataL,
                             pdataR,
                             dir=None):  # centered flux ; pL[ieq][face]
        gam = self.gamma
        gam1 = gam - 1.

        rhoL = pdataL[0]
        uL = unL = pdataL[1]
        pL = pdataL[2]
        rhoR = pdataR[0]
        uR = unR = pdataR[1]
        pR = pdataR[2]

        cL2 = gam * pL / rhoL
        cR2 = gam * pR / rhoR
        HL = cL2 / gam1 + 0.5 * uL**2
        HR = cR2 / gam1 + 0.5 * uR**2

        # final flux
        Frho = .5 * (rhoL * unL + rhoR * unR)
        Frhou = .5 * ((rhoL * unL**2 + pL) + (rhoR * unR**2 + pR))
        FrhoE = .5 * ((rhoL * unL * HL) + (rhoR * unR * HR))

        return [Frho, Frhou, FrhoE]

    @_numfluxdict.register()
    def numflux_centeredmassflow(self,
                                 pdataL,
                                 pdataR,
                                 dir=None):  # centered flux ; pL[ieq][face]
        gam = self.gamma
        gam1 = gam - 1.

        rhoL = pdataL[0]
        uL = unL = pdataL[1]
        pL = pdataL[2]
        rhoR = pdataR[0]
        uR = unR = pdataR[1]
        pR = pdataR[2]

        cL2 = gam * pL / rhoL
        cR2 = gam * pR / rhoR
        HL = cL2 / gam1 + 0.5 * uL**2
        HR = cR2 / gam1 + 0.5 * uR**2

        # final flux
        Frho = .5 * (rhoL * unL + rhoR * unR)
        Frhou = .5 * (Frho * (unL + unR) + pL + pR)
        FrhoE = .5 * Frho * (HL + HR)

        return [Frho, Frhou, FrhoE]

    @_numfluxdict.register()
    def numflux_hlle(self,
                     pdataL,
                     pdataR,
                     dir=None):  # HLLE Riemann solver ; pL[ieq][face]

        gam = self.gamma
        gam1 = gam - 1.

        rhoL = pdataL[0]
        uL = unL = pdataL[1]
        pL = pdataL[2]
        rhoR = pdataR[0]
        uR = unR = pdataR[1]
        pR = pdataR[2]

        cL2 = gam * pL / rhoL
        cR2 = gam * pR / rhoR
        HL = cL2 / gam1 + 0.5 * uL**2
        HR = cR2 / gam1 + 0.5 * uR**2

        # The HLLE Riemann solver

        # sorry for using little "e" here - is is not just internal energy
        eL = HL - pL / rhoL
        eR = HR - pR / rhoR

        # Roe's averaging
        Rrho, uRoe, cRoe = self._Roe_average(rhoL, uL, HL, rhoR, uR, HR)
        # max HLL 2 waves "velocities"
        sL = np.minimum(0., np.minimum(uRoe - cRoe, unL - np.sqrt(cL2)))
        sR = np.maximum(0., np.maximum(uRoe + cRoe, unR + np.sqrt(cR2)))

        # final flux
        Frho = (sR * rhoL * unL - sL * rhoR * unR + sL * sR *
                (rhoR - rhoL)) / (sR - sL)
        Frhou = (sR * (rhoL * unL**2 + pL) - sL *
                 (rhoR * unR**2 + pR) + sL * sR *
                 (rhoR * unR - rhoL * unL)) / (sR - sL)
        FrhoE = (sR * (rhoL * unL * HL) - sL * (rhoR * unR * HR) + sL * sR *
                 (rhoR * eR - rhoL * eL)) / (sR - sL)

        return [Frho, Frhou, FrhoE]

    @_numfluxdict.register()
    def numflux_hllc(self,
                     pdataL,
                     pdataR,
                     dir=None):  # HLLC Riemann solver ; pL[ieq][face]

        gam = self.gamma
        gam1 = gam - 1.

        rhoL = pdataL[0]
        uL = unL = pdataL[1]
        pL = pdataL[2]
        rhoR = pdataR[0]
        uR = unR = pdataR[1]
        pR = pdataR[2]

        cL2 = gam * pL / rhoL
        cR2 = gam * pR / rhoR

        # the enthalpy is assumed to include ke ...!
        HL = cL2 / gam1 + 0.5 * uL**2
        HR = cR2 / gam1 + 0.5 * uR**2

        # The HLLC Riemann solver

        # sorry for using little "e" here - is is not just internal energy
        eL = HL - pL / rhoL
        eR = HR - pR / rhoR

        # Roe's averaging
        Rrho, uRoe, cRoe = self._Roe_average(rhoL, uL, HL, rhoR, uR, HR)

        # speed of sound at L and R
        sL = np.minimum(uRoe - cRoe, unL - np.sqrt(cL2))
        sR = np.maximum(uRoe + cRoe, unR + np.sqrt(cR2))

        # speed of contact surface
        sM = (pL - pR - rhoL * unL * (sL - unL) + rhoR * unR *
              (sR - unR)) / (rhoR * (sR - unR) - rhoL * (sL - unL))

        # pressure at right and left (pR=pL) side of contact surface
        pStar = rhoR * (unR - sR) * (unR - sM) + pR

        # should not be computed if totally upwind
        SmoSSm = np.where(sM >= 0., sM / (sL - sM), sM / (sR - sM))
        SmUoSSm = np.where(sM >= 0., (sL - unL) / (sL - sM),
                           (sR - unR) / (sR - sM))

        Frho = np.where(sM >= 0.,
                        np.where(sL >= 0., rhoL * unL, rhoL * sM * SmUoSSm),
                        np.where(sR <= 0., rhoR * unR, rhoR * sM * SmUoSSm))

        Frhou = np.where(
            sM >= 0.,
            np.where(sL >= 0., Frho * uL + pL,
                     Frho * uL + (pStar - pL) * SmoSSm + pStar),
            np.where(sR <= 0., Frho * uR + pR,
                     Frho * uR + (pStar - pR) * SmoSSm + pStar))

        FrhoE = np.where(
            sM >= 0.,
            np.where(sL >= 0., rhoL * HL * unL, Frho * eL +
                     (pStar * sM - pL * unL) * SmoSSm + pStar * sM),
            np.where(sR <= 0., rhoR * HR * unR, Frho * eR +
                     (pStar * sM - pR * unR) * SmoSSm + pStar * sM))

        return [Frho, Frhou, FrhoE]

    def timestep(self, data, dx, condition):
        "computation of timestep with conservative data"
        #        dt = CFL * dx / ( |u| + c )
        # dt = np.zeros(len(dx)) #test use zeros instead
        #dt = condition*dx/ (data[1] + np.sqrt(self.gamma*data[2]/data[0]) )
        Vmag = self.velocitymag(data)
        dt = condition * dx / (Vmag +
                               np.sqrt(self.gamma * (self.gamma - 1.0) *
                                       (data[2] / data[0] - 0.5 * Vmag**2)))
        return dt
Exemple #9
0
class shallowwater1d(base.model):
    """
    Class model for shallow water equations

    attributes:

    """
    _bcdict = base.methoddict(
        'bc_')  # dict and associated decorator method to register BC
    _vardict = base.methoddict()
    _numfluxdict = base.methoddict('numflux_')

    def __init__(self, g=9.81, source=None):
        base.model.__init__(self, name='shallowwater', neq=2)
        self.islinear = 0
        self.shape = [1, 1]
        self.g = g  # gravity attraction
        self.source = source
        self._bcdict.merge(shallowwater1d._bcdict)
        self._vardict.merge(shallowwater1d._vardict)
        self._numfluxdict.merge(shallowwater1d._numfluxdict)

    def cons2prim(self, qdata):  # qdata[ieq][cell] :
        """
        >>> model().cons2prim([[5.], [10.], [20.]])
        True
        """
        h = qdata[0]
        u = qdata[1] / qdata[0]
        pdata = [h, u]
        return pdata  #p pdata : primitive variables

    def prim2cons(self, pdata):  # qdata[ieq][cell] :
        """
        >>> model().prim2cons([[2.], [4.], [10.]]) == [[2.], [8.], [41.]]
        True
        """
        qdata = [pdata[0], pdata[0] * pdata[1]]
        return qdata  # Convervative variables

    @_vardict.register()
    def height(self, qdata):
        return qdata[0].copy()

    @_vardict.register()
    def massflow(self, qdata):
        return qdata[1].copy()

    @_vardict.register()
    def velocity(self, qdata):
        return qdata[1] / qdata[0]

    def numflux(self, name, pdataL, pdataR, dir=None):
        if name is None: name = 'rusanov'
        return (self._numfluxdict.dict[name])(self, pdataL, pdataR, dir)

    @_numfluxdict.register(name='centered')
    @_numfluxdict.register()
    def numflux_centeredflux(
        self,
        pdataL,
        pdataR,
        dir=None
    ):  # centered flux ; pL[ieq][face] : in primitive variables (h,u) !
        g = self.g

        hL = pdataL[0]
        uL = pdataL[1]
        hR = pdataR[0]
        uR = pdataR[1]

        # final flux
        Fh = .5 * (hL * uL + hR * uR)
        Fq = .5 * ((hL * uL**2 + 0.5 * g * hL**2) +
                   (hR * uR**2 + 0.5 * g * hR**2))

        return [Fh, Fq]

    @_numfluxdict.register()
    def numflux_rusanov(self,
                        pdataL,
                        pdataR,
                        dir=None):  # Rusanov flux ; pL[ieq][face]
        g = self.g
        #
        hL = pdataL[0]
        uL = pdataL[1]
        hR = pdataR[0]
        uR = pdataR[1]
        cL = np.sqrt(g * hL)
        cR = np.sqrt(g * hR)
        # Eigen values definition
        cmax = np.maximum(abs(uL) + cL, abs(uR) + cR)
        # final Rusanov flux
        qL = hL * uL
        qR = hR * uR
        Fh = .5 * (qL + qR) - 0.5 * cmax * (hR - hL)
        Fq = .5 * ((qL * uL + 0.5 * g * hL**2) +
                   (qR * uR**2 + 0.5 * g * hR**2)) - 0.5 * cmax * (qR - qL)

        return [Fh, Fq]

    @_numfluxdict.register()
    def numflux_hll(self,
                    pdataL,
                    pdataR,
                    dir=None):  # HLL flux ; pL[ieq][face]
        g = self.g

        hL = pdataL[0]
        uL = pdataL[1]
        hR = pdataR[0]
        uR = pdataR[1]

        cL = np.sqrt(g * hL)
        cR = np.sqrt(g * hR)
        sL = np.minimum(0., np.minimum(uL - cL, uR - cR))
        sR = np.maximum(0., np.maximum(uL + cL, uR + cR))
        kL = sR / (sR - sL)
        kR = -sL / (sR - sL)
        qL = hL * uL
        qR = hR * uR
        Fh = kL * qL + kR * qR - kR * sR * (hR - hL)
        Fu = kL * (qL * uL + .5 * g * hL**2) + kR * (
            qR * uR + .5 * g * hR**2) - kR * sR * (qR - qL)
        return [Fh, Fu]

    def timestep(self, data, dx, condition):
        "computation of timestep: data(=pdata) is not used, dx is an array of cell sizes, condition is the CFL number"
        #        dt = CFL * dx / c where c is the highest eigen value velocity
        g = self.g
        c = abs(data[1] / data[0]) + np.sqrt(g * data[0])
        dt = condition * dx / c
        return dt

    @_bcdict.register()
    def bc_sym(self, dir, data, param):  # In primitive values here
        "symmetry boundary condition, for inviscid equations, it is equivalent to a wall, do not need user parameters"
        return [data[0], -data[1]]

    @_bcdict.register()
    def bc_inf(self, dir, data,
               param):  # Simulate infinite plane : not sure...
        #zeros_h = np.zeros( np.shape(data[0]))
        #zeros_u = np.zeros( np.shape(data[1]))
        return [data[0], data[1]]