Пример #1
0
 def post_loop(self, d_idx, d_m, d_u, d_v, d_h, d_closest_idx, d_merge, d_x,
               d_y, KERNEL):
     idx = declare('int')
     xma = declare('matrix(3)')
     xmb = declare('matrix(3)')
     idx = d_closest_idx[d_idx]
     if idx > -1:
         if d_idx == d_closest_idx[idx]:
             if d_idx < idx:
                 m_merged = d_m[d_idx] + d_m[idx]
                 x_merged = (d_m[d_idx]*d_x[d_idx] + d_m[idx]*d_x[idx]) \
                             / m_merged
                 y_merged = (d_m[d_idx]*d_y[d_idx] + d_m[idx]*d_y[idx]) \
                             / m_merged
                 xma[0] = x_merged - d_x[d_idx]
                 xma[1] = y_merged - d_y[d_idx]
                 xmb[0] = x_merged - d_x[idx]
                 xmb[1] = y_merged - d_y[idx]
                 rma = sqrt(xma[0] * xma[0] + xma[1] * xma[1])
                 rmb = sqrt(xmb[0] * xmb[0] + xmb[1] * xmb[1])
                 d_u[d_idx] = (d_m[d_idx]*d_u[d_idx] + d_m[idx]*d_u[idx]) \
                               / m_merged
                 d_v[d_idx] = (d_m[d_idx]*d_v[d_idx] + d_m[idx]*d_v[idx]) \
                               / m_merged
                 const1 = d_m[d_idx] * KERNEL.kernel(xma, rma, d_h[d_idx])
                 const2 = d_m[idx] * KERNEL.kernel(xmb, rmb, d_h[idx])
                 d_h[d_idx] = sqrt(
                     (7 * pi / 10.) * (m_merged / (const1 + const2)))
                 d_m[d_idx] = m_merged
             else:
                 d_merge[d_idx] = 1
Пример #2
0
def gj_solve(A=[1., 0.], b=[1., 0.], n=3, result=[0., 1.]):
    r""" A gauss-jordan method to solve an augmented matrix for the
        unknown variables, x, in Ax = b.

    References
    ----------
    https://ricardianambivalence.com/2012/10/20/pure-python-gauss-jordan
    -solve-ax-b-invert-a/
    """

    i, j, eqns, colrange, augCol, col, row, bigrow = declare('int', 8)
    m = declare('matrix((4, 5))')
    for i in range(n):
        for j in range(n):
            m[i][j] = A[n * i + j]
        m[i][n] = b[i]

    eqns = n
    colrange = n
    augCol = n + 1

    for col in range(colrange):
        bigrow = col
        for row in range(col + 1, colrange):
            if abs(m[row][col]) > abs(m[bigrow][col]):
                bigrow = row
                temp = m[row][col]
                m[row][col] = m[bigrow][col]
                m[bigrow][col] = temp

    rr, rrcol, rb, rbr, kup, kupr, kleft, kleftr = declare('int', 8)
    for rrcol in range(0, colrange):
        for rr in range(rrcol + 1, eqns):
            cc = -(float(m[rr][rrcol]) / float(m[rrcol][rrcol]))
            for j in range(augCol):
                m[rr][j] = m[rr][j] + cc * m[rrcol][j]

    backCol, backColr = declare('int', 2)
    tol = 1.0e-05
    for rbr in range(eqns):
        rb = eqns - rbr - 1
        if (m[rb][rb] == 0):
            if abs(m[rb][augCol - 1]) >= tol:
                return 0.0
        else:
            for backColr in range(rb, augCol):
                backCol = rb + augCol - backColr - 1
                m[rb][backCol] = float(m[rb][backCol]) / float(m[rb][rb])
            if not (rb == 0):
                for kupr in range(rb):
                    kup = rb - kupr - 1
                    for kleftr in range(rb, augCol):
                        kleft = rb + augCol - kleftr - 1
                        kk = -float(m[kup][rb]) / float(m[rb][rb])
                        m[kup][kleft] = m[kup][kleft] + \
                            kk * float(m[rb][kleft])
    for i in range(n):
        result[i] = m[i][n]
Пример #3
0
 def loop_all(self, d_sum_Ak, d_pa_to_split, d_parent_idx, d_idx, s_m,
              s_rho, s_parent_idx, NBRS, N_NBRS):
     i = declare('int')
     s_idx = declare('long')
     if d_pa_to_split[d_idx]:
         for i in range(N_NBRS):
             s_idx = NBRS[i]
             if s_parent_idx[s_idx] == d_parent_idx[d_idx]:
                 d_sum_Ak[d_idx] += s_m[s_idx] / s_rho[s_idx]
Пример #4
0
def roe(rhol=0.0,
        rhor=1.0,
        pl=0.0,
        pr=1.0,
        ul=0.0,
        ur=1.0,
        gamma=1.4,
        niter=20,
        tol=1e-6,
        result=[0.0, 0.0]):
    """Roe's approximate Riemann solver for the Euler equations to
    determine the intermediate states.

    Parameters
    ----------
    rhol, rhor: double: left and right density.
    pl, pr: double: left and right pressure.
    ul, ur: double: left and right speed.
    gamma: double: Ratio of specific heats.
    niter: int: Max number of iterations to try for iterative schemes.
    tol: double: Error tolerance for convergence.

    result: list: List of length 2. Result will contain computed pstar, ustar

    Returns
    -------

    Returns 0 if the value is computed correctly else 1 if the iterations
    either did not converge or if there is an error.

    """

    rrhol, rrhor, denominator, plr, vlr, ulr = declare('double', 6)
    cslr, cslr1 = declare('double', 2)
    # Roe-averaged pressure and specific volume
    rrhol = sqrt(rhol)
    rrhor = sqrt(rhor)

    denominator = 1. / (rrhor + rrhol)

    plr = (rrhol * pl + rrhor * pr) * denominator
    vlr = (rrhol / rhol + rrhor / rhor) * denominator
    ulr = (rrhol * ul + rrhor * ur) * denominator

    # average sound speed at the interface
    cslr = sqrt(gamma * plr / vlr)
    cslr1 = 1. / cslr

    # the intermediate states
    result[0] = plr - 0.5 * (ur - ul) * cslr
    result[1] = ulr - 0.5 * (pr - pl) * cslr1

    return 0
Пример #5
0
 def loop_all(self, d_shep_corr, d_x, d_y, d_idx, s_x, s_y, s_m, s_rho,
              s_idx, d_h, KERNEL, NBRS, N_NBRS):
     i = declare('int')
     xij = declare('matrix(3)')
     rij = 0.0
     corr_sum = 0.0
     for i in range(N_NBRS):
         s_idx = NBRS[i]
         xij[0] = d_x[d_idx] - s_x[s_idx]
         xij[1] = d_y[d_idx] - s_y[s_idx]
         rij = sqrt(xij[0] * xij[0] + xij[1] * xij[1])
         corr_sum += (s_m[s_idx] / s_rho[s_idx]) * KERNEL.kernel(
             xij, rij, d_h[d_idx])
     d_shep_corr[d_idx] = corr_sum
Пример #6
0
 def loop_all(self, d_idx, d_x, d_rho, s_m, s_x, s_h, KERNEL, NBRS, N_NBRS):
     i = declare('int')
     s_idx = declare('long')
     xij = declare('matrix((3,))')
     rij = 0.0
     sum = 0.0
     xij[1] = 0.0
     xij[2] = 0.0
     for i in range(N_NBRS):
         s_idx = NBRS[i]
         xij[0] = d_x[d_idx] - s_x[s_idx]
         rij = fabs(xij[0])
         sum += s_m[s_idx]*KERNEL.kernel(xij, rij, s_h[s_idx])
     d_rho[d_idx] += sum
Пример #7
0
 def loop_all(self, d_idx, d_rho, d_x, d_y, d_z, s_x, s_y, s_z, d_h, s_h,
              s_m, s_rhotmp, KERNEL, NBRS, N_NBRS):
     n, i, j, k, s_idx = declare('int', 5)
     n = 4
     amls = declare('matrix(16)')
     x = d_x[d_idx]
     y = d_y[d_idx]
     z = d_z[d_idx]
     xij = declare('matrix(4)')
     for i in range(n):
         for j in range(n):
             amls[n * i + j] = 0.0
     for k in range(N_NBRS):
         s_idx = NBRS[k]
         xij[0] = x - s_x[s_idx]
         xij[1] = y - s_y[s_idx]
         xij[2] = z - s_z[s_idx]
         rij = sqrt(xij[0] * xij[0] + xij[1] * xij[1] + xij[2] * xij[2])
         hij = (d_h[d_idx] + s_h[s_idx]) * 0.5
         wij = KERNEL.kernel(xij, rij, hij)
         for i in range(n):
             if i == 0:
                 fac1 = 1.0
             else:
                 fac1 = xij[i - 1]
             for j in range(n):
                 if j == 0:
                     fac2 = 1.0
                 else:
                     fac2 = xij[j - 1]
                 amls[n * i + j] += fac1 * fac2 * \
                     s_m[s_idx] * wij / s_rhotmp[s_idx]
     res = declare('matrix(4)')
     gj_solve(amls, [1., 0., 0., 0.], n, res)
     b0 = res[0]
     b1 = res[1]
     b2 = res[2]
     b3 = res[3]
     d_rho[d_idx] = 0.0
     for k in range(N_NBRS):
         s_idx = NBRS[k]
         xij[0] = x - s_x[s_idx]
         xij[1] = y - s_y[s_idx]
         xij[2] = z - s_z[s_idx]
         rij = sqrt(xij[0] * xij[0] + xij[1] * xij[1] + xij[2] * xij[2])
         hij = (d_h[d_idx] + s_h[s_idx]) * 0.5
         wij = KERNEL.kernel(xij, rij, hij)
         wmls = (b0 + b1 * xij[0] + b2 * xij[1] + b3 * xij[2]) * wij
         d_rho[d_idx] += s_m[s_idx] * wmls
Пример #8
0
 def loop_all(self, d_rho, d_idx, d_merge, d_closest_idx, d_x, d_y, s_h,
              s_m, s_x, s_y, KERNEL, NBRS, N_NBRS):
     i = declare('int')
     s_idx = declare('long')
     xij = declare('matrix(3)')
     if d_merge[d_closest_idx[d_idx]] == 1:
         d_rho[d_idx] = 0.0
         rij = 0.0
         rho_sum = 0.0
         for i in range(N_NBRS):
             s_idx = NBRS[i]
             xij[0] = d_x[d_idx] - s_x[s_idx]
             xij[1] = d_y[d_idx] - s_y[s_idx]
             rij = sqrt(xij[0] * xij[0] + xij[1] * xij[1])
             rho_sum += s_m[s_idx] * KERNEL.kernel(xij, rij, s_h[s_idx])
         d_rho[d_idx] += rho_sum
Пример #9
0
def prefun_exact(p=0.0,
                 dk=0.0,
                 pk=0.0,
                 ck=0.0,
                 g1=0.0,
                 g2=0.0,
                 g4=0.0,
                 g5=0.0,
                 g6=0.0,
                 result=[0.0, 0.0]):
    """The pressure function.  Updates result with f, fd.
    """

    pratio, f, fd, ak, bk, qrt = declare('double', 6)

    if (p <= pk):
        pratio = p / pk
        f = g4 * ck * (pratio**g1 - 1.0)
        fd = (1.0 / (dk * ck)) * pratio**(-g2)
    else:
        ak = g5 / dk
        bk = g6 * pk
        qrt = sqrt(ak / (bk + p))
        f = (p - pk) * qrt
        fd = (1.0 - 0.5 * (p - pk) / (bk + p)) * qrt

    result[0] = f
    result[1] = fd
Пример #10
0
    def reduce(self, dst):
        n = len(dst.x)
        tmp_sum_logrho = serial_reduce_array(dst.logrho, 'sum')
        sum_logrho = parallel_reduce_array(tmp_sum_logrho, 'sum')
        g = exp(sum_logrho / n)

        lamda = declare('object')
        lamda = self.k * numpy.power(g / dst.rho, self.eps)
        dst.h[:] = lamda * dst.h0
Пример #11
0
 def loop_all(self, d_idx, d_rho, d_x, d_y, d_z, s_m, s_rhotmp, s_x, s_y,
              s_z, d_h, s_h, KERNEL, NBRS, N_NBRS):
     i, s_idx = declare('int', 2)
     xij = declare('matrix(3)')
     tmp_w = 0.0
     x = d_x[d_idx]
     y = d_y[d_idx]
     z = d_z[d_idx]
     d_rho[d_idx] = 0.0
     for i in range(N_NBRS):
         s_idx = NBRS[i]
         xij[0] = x - s_x[s_idx]
         xij[1] = y - s_y[s_idx]
         xij[2] = z - s_z[s_idx]
         rij = sqrt(xij[0] * xij[0] + xij[1] * xij[1] + xij[2] * xij[2])
         hij = (d_h[d_idx] + s_h[s_idx]) * 0.5
         wij = KERNEL.kernel(xij, rij, hij)
         tmp_w += wij * s_m[s_idx] / s_rhotmp[s_idx]
         d_rho[d_idx] += wij * s_m[s_idx]
     d_rho[d_idx] /= tmp_w
Пример #12
0
def llxf(rhol=0.0,
         rhor=1.0,
         pl=0.0,
         pr=1.0,
         ul=0.0,
         ur=1.0,
         gamma=1.4,
         niter=20,
         tol=1e-6,
         result=[0.0, 0.0]):
    """Lax Friedrichs approximate Riemann solver for the Euler equations to
    determine the intermediate states.

    Parameters
    ----------
    rhol, rhor: double: left and right density.
    pl, pr: double: left and right pressure.
    ul, ur: double: left and right speed.
    gamma: double: Ratio of specific heats.
    niter: int: Max number of iterations to try for iterative schemes.
    tol: double: Error tolerance for convergence.

    result: list: List of length 2. Result will contain computed pstar, ustar

    Returns
    -------

    Returns 0 if the value is computed correctly else 1 if the iterations
    either did not converge or if there is an error.

    """

    gamma1, csl, csr, cslr, el, El, er, Er, pstar = declare('double', 9)
    gamma1 = 1. / (gamma - 1.0)

    # Lagrangian sound speeds
    csl = sqrt(gamma * pl * rhol)
    csr = sqrt(gamma * pr * rhor)
    cslr = max(csr, csl)

    # Total energy on either side
    el = pl * gamma1 / rhol
    El = el + 0.5 * ul * ul

    er = pr * gamma1 / rhor
    Er = er + 0.5 * ur * ur

    # the intermediate states
    # cdef double ustar = 0.5 * ( ul + ur - 1./cslr * (pr - pl) )
    pstar = 0.5 * (pl + pr - cslr * (ur - ul))
    result[0] = pstar
    result[1] = 1. / pstar * (0.5 * ((pl * ul + pr * ur) - cslr * (Er - El)))
    return 0
Пример #13
0
 def loop_all(self, d_idx, d_merge, d_closest_idx, d_x, d_y, d_h, d_A, s_x,
              s_y, s_A, NBRS, N_NBRS):
     i, closest = declare('int', 2)
     s_idx = declare('unsigned int')
     d_merge[d_idx] = 0
     xi = d_x[d_idx]
     yi = d_y[d_idx]
     rmin = d_h[d_idx] * 10.0
     closest = -1
     if (d_A[d_idx] < self.A_min
             and ((self.x_min < d_x[d_idx] < self.x_max) and
                  (self.y_min < d_y[d_idx] < self.y_max))):
         for i in range(N_NBRS):
             s_idx = NBRS[i]
             if s_idx == d_idx:
                 continue
             xij = xi - s_x[s_idx]
             yij = yi - s_y[s_idx]
             rij = sqrt(xij * xij + yij * yij)
             if rij < rmin:
                 closest = s_idx
                 rmin = rij
     d_closest_idx[d_idx] = closest
Пример #14
0
def hllc(rhol=0.0,
         rhor=1.0,
         pl=0.0,
         pr=1.0,
         ul=0.0,
         ur=1.0,
         gamma=1.4,
         niter=20,
         tol=1e-6,
         result=[0.0, 0.0]):
    """Harten-Lax-van Leer-Contact approximate Riemann solver for the Euler
    equations to determine the intermediate states.

    Parameters
    ----------
    rhol, rhor: double: left and right density.
    pl, pr: double: left and right pressure.
    ul, ur: double: left and right speed.
    gamma: double: Ratio of specific heats.
    niter: int: Max number of iterations to try for iterative schemes.
    tol: double: Error tolerance for convergence.

    result: list: List of length 2. Result will contain computed pstar, ustar

    Returns
    -------

    Returns 0 if the value is computed correctly else 1 if the iterations
    either did not converge or if there is an error.

    """

    gamma1, pstar, ustar, mstar, estar = declare('double', 5)
    gamma1 = 1. / (gamma - 1.0)

    # Roe-averaged interface velocity
    rrhol = sqrt(rhol)
    rrhor = sqrt(rhor)
    ulr = (rrhol * ul + rrhor * ur) / (rrhol + rrhor)

    # relative velocities at the interface
    vl = ul - ulr
    vr = ur - ulr

    # Sound speeds at the interface
    csl = sqrt(gamma * pl / rhol)
    csr = sqrt(gamma * pr / rhor)
    cslr = (rrhol * csl + rrhor * csr) / (rrhol + rrhor)

    # wave speed estimates at the interface
    sl = min(vl - csl, 0 - cslr)
    sr = max(vr + csr, 0 + cslr)

    sm = rhor * vr * (sr - vr) - rhol * vl * (sl - vl) + pl - pr
    sm /= (rhor * (sr - vr) - rhol * (sl - vl))

    # phat
    phat = rhol * (vl - sl) * (vl - sm) + pl

    # Total energy on either side
    el = pl * gamma1 / rhol
    El = rhol * (el + 0.5 * ul * ul)

    er = pr * gamma1 / rhor
    Er = rhor * (er + 0.5 * ur * ur)

    # Momentum on either side
    Ml = rhol * ul
    Mr = rhor * ur

    # compute the values based on the wave speeds
    if (sl > 0):
        pstar = pl
        ustar = ul

    elif (sl <= 0.0) and (0.0 < sm):
        mstar = 1. / (sl - sm) * ((sl - vl) * Ml + (phat - pl))
        estar = 1. / (sl - sm) * ((sl - vl) * El - pl * vl + phat * sm)

        pstar = sm * mstar + phat

        ustar = sm * estar + (sm + ulr) * phat
        ustar /= pstar

    elif (sm <= 0) and (0 < sr):
        mstar = 1. / (sr - sm) * ((sr - vr) * Mr + (phat - pr))
        estar = 1. / (sr - sm) * ((sr - vr) * Er - pr * vr + phat * sm)

        pstar = sm * mstar + phat

        ustar = sm * estar + (sm + ulr) * phat
        ustar /= pstar

    elif (sr < 0):
        pstar = pr
        ustar = ur

    else:
        printf("Incorrect wave speeds")
        return 1

    result[0] = pstar
    result[1] = ustar
    return 0
Пример #15
0
 def func(self, d_idx, d_x=[0.0, 0.0]):
     val, val1 = declare('float', 2)
     # val1 = declare('double')
     val = d_x[d_idx]
     index = declare('unsigned int')
     index = d_idx
Пример #16
0
def exact(rhol=0.0,
          rhor=1.0,
          pl=0.0,
          pr=1.0,
          ul=0.0,
          ur=1.0,
          gamma=1.4,
          niter=20,
          tol=1e-6,
          result=[0.0, 0.0]):
    """Exact Riemann solver.

    Solve the Riemann problem for the Euler equations to determine the
    intermediate states.

    Parameters
    ----------
    rhol, rhor: double: left and right density.
    pl, pr: double: left and right pressure.
    ul, ur: double: left and right speed.
    gamma: double: Ratio of specific heats.
    niter: int: Max number of iterations to try for iterative schemes.
    tol: double: Error tolerance for convergence.

    result: list: List of length 2. Result will contain computed pstar, ustar

    Returns
    -------

    Returns 0 if the value is computed correctly else 1 if the iterations
    either did not converge or if there is an error.

    """
    # save derived variables
    tmp1 = 1.0 / (2 * gamma)
    tmp2 = 1.0 / (gamma - 1.0)
    tmp3 = 1.0 / (gamma + 1.0)

    gamma1 = (gamma - 1.0) * tmp1
    gamma2 = (gamma + 1.0) * tmp1
    gamma3 = 2 * gamma * tmp2
    gamma4 = 2 * tmp2
    gamma5 = 2 * tmp3
    gamma6 = tmp3 / tmp2
    gamma7 = 0.5 * (gamma - 1.0)

    # sound speed to the left and right based on the Ideal EOS
    cl, cr = declare('double', 2)
    cl = sqrt(gamma * pl / rhol)
    cr = sqrt(gamma * pr / rhor)

    # Iteration constants
    i = declare('int')
    i = 0

    # local variables
    qser, cup, ppv, pmin, pmax, qmax = declare('double', 6)
    pq, um, ptl, ptr, pm, gel, ger = declare('double', 7)
    pstar, pold, udifff = declare('double', 3)
    p, change = declare('double', 2)

    # initialize variables
    fl, fr = declare('matrix(2)', 2)
    p = 0.0

    # check the initial data
    if (gamma4 * (cl + cr) <= (ur - ul)):
        return 1

    # compute the guess pressure 'pm'
    qser = 2.0
    cup = 0.25 * (rhol + rhor) * (cl + cr)
    ppv = 0.5 * (pl + pr) + 0.5 * (ul - ur) * cup
    pmin = min(pl, pr)
    pmax = max(pl, pr)
    qmax = pmax / pmin

    ppv = max(0.0, ppv)

    if ((qmax <= qser) and (pmin <= ppv) and (ppv <= pmax)):
        pm = ppv
    elif (ppv < pmin):
        pq = (pl / pr)**gamma1
        um = (pq * ul / cl + ur / cr + gamma4 *
              (pq - 1.0)) / (pq / cl + 1.0 / cr)
        ptl = 1.0 + gamma7 * (ul - um) / cl
        ptr = 1.0 + gamma7 * (um - ur) / cr
        pm = 0.5 * (pl * ptl**gamma3 + pr * ptr**gamma3)
    else:
        gel = sqrt((gamma5 / rhol) / (gamma6 * pl + ppv))
        ger = sqrt((gamma5 / rhor) / (gamma6 * pr + ppv))
        pm = (gel * pl + ger * pr - (ur - ul)) / (gel + ger)

    # the guessed value is pm
    pstart = pm

    pold = pstart
    udifff = ur - ul

    for i in range(niter):
        prefun_exact(pold, rhol, pl, cl, gamma1, gamma2, gamma4, gamma5,
                     gamma6, fl)
        prefun_exact(pold, rhor, pr, cr, gamma1, gamma2, gamma4, gamma5,
                     gamma6, fr)

        p = pold - (fl[0] + fr[0] + udifff) / (fl[1] + fr[1])
        change = 2.0 * abs((p - pold) / (p + pold))
        if change <= tol:
            break
        pold = p

    if i == niter - 1:
        printf("Divergence in Newton-Raphson Iteration")
        return 1

    # compute the velocity in the star region 'um'
    um = 0.5 * (ul + ur + fr[0] - fl[0])
    result[0] = p
    result[1] = um
    return 0
Пример #17
0
def sample(pm=0.0,
           um=0.0,
           s=0.0,
           rhol=1.0,
           rhor=0.0,
           pl=1.0,
           pr=0.0,
           ul=1.0,
           ur=0.0,
           gamma=1.4,
           result=[0.0, 0.0, 0.0]):
    """Sample the solution to the Riemann problem.

    Parameters
    ----------

    pm, um : float
        Pressure and velocity in the star region as returned
        by `exact`

    s : float
        Sampling point (x/t)

    rhol, rhor : float
        Densities on either side of the discontinuity

    pl, pr : float
        Pressures on either side of the discontinuity

    ul, ur : float
        Velocities on either side of the discontinuity

    gamma : float
        Ratio of specific heats

    result : list(double)
        (rho, u, p) : Sampled density, velocity and pressure

    """
    tmp1, tmp2, tmp3 = declare('double', 3)
    g1, g2, g3, g4, g5, g6 = declare('double', 6)
    cl, cr = declare('double', 2)

    # save derived variables
    tmp1 = 1.0 / (2 * gamma)
    tmp2 = 1.0 / (gamma - 1.0)
    tmp3 = 1.0 / (gamma + 1.0)

    g1 = (gamma - 1.0) * tmp1
    g2 = (gamma + 1.0) * tmp1
    g3 = 2 * gamma * tmp2
    g4 = 2 * tmp2
    g5 = 2 * tmp3
    g6 = tmp3 / tmp2
    g7 = 0.5 * (gamma - 1.0)

    # sound speeds at the left and right data states
    cl = sqrt(gamma * pl / rhol)
    cr = sqrt(gamma * pr / rhor)

    if s <= um:
        # sampling point lies to the left of the contact discontinuity
        if (pm <= pl):  # left rarefaction
            # speed of the head of the rarefaction
            shl = ul - cl

            if (s <= shl):
                # sampled point is left state
                rho = rhol
                u = ul
                p = pl
            else:
                cml = cl * (pm / pl)**g1  # eqn (4.54)
                stl = um - cml  # eqn (4.55)

                if (s > stl):
                    # sampled point is Star left state. eqn (4.53)
                    rho = rhol * (pm / pl)**(1.0 / gamma)
                    u = um
                    p = pm
                else:
                    # sampled point is inside left fan
                    u = g5 * (cl + g7 * ul + s)
                    c = g5 * (cl + g7 * (ul - s))

                    rho = rhol * (c / cl)**g4
                    p = pl * (c / cl)**g3

        else:  # pm <= pl
            # left shock
            pml = pm / pl
            sl = ul - cl * sqrt(g2 * pml + g1)

            if (s <= sl):
                # sampled point is left data state
                rho = rhol
                u = ul
                p = pl
            else:
                # sampled point is Star left state
                rho = rhol * (pml + g6) / (pml * g6 + 1.0)
                u = um
                p = pm

    else:  # s < um
        # sampling point lies to the right of the contact discontinuity
        if (pm > pr):
            # right shock
            pmr = pm / pr
            sr = ur + cr * sqrt(g2 * pmr + g1)

            if (s >= sr):
                # sampled point is right data state
                rho = rhor
                u = ur
                p = pr
            else:
                # sampled point is star right state
                rho = rhor * (pmr + g6) / (pmr * g6 + 1.0)
                u = um
                p = pm
        else:
            # right rarefaction
            shr = ur + cr

            if (s >= shr):
                # sampled point is right state
                rho = rhor
                u = ur
                p = pr
            else:
                cmr = cr * (pm / pr)**g1
                STR = um + cmr

                if (s <= STR):
                    # sampled point is star right
                    rho = rhor * (pm / pr)**(1.0 / gamma)
                    u = um
                    p = pm
                else:
                    # sampled point is inside left fan
                    u = g5 * (-cr + g7 * ur + s)
                    c = g5 * (cr - g7 * (ur - s))
                    rho = rhor * (c / cr)**g4
                    p = pr * (c / cr)**g3

    result[0] = rho
    result[1] = u
    result[2] = p
Пример #18
0
def ducowicz(rhol=0.0,
             rhor=1.0,
             pl=0.0,
             pr=1.0,
             ul=0.0,
             ur=1.0,
             gamma=1.4,
             niter=20,
             tol=1e-6,
             result=[0.0, 0.0]):
    """Ducowicz Approximate Riemann solver for the Euler equations to
    determine the intermediate states.

    Parameters
    ----------
    rhol, rhor: double: left and right density.
    pl, pr: double: left and right pressure.
    ul, ur: double: left and right speed.
    gamma: double: Ratio of specific heats.
    niter: int: Max number of iterations to try for iterative schemes.
    tol: double: Error tolerance for convergence.

    result: list: List of length 2. Result will contain computed pstar, ustar

    Returns
    -------

    Returns 0 if the value is computed correctly else 1 if the iterations
    either did not converge or if there is an error.

    """
    al, ar = declare('double', 2)
    # strong shock parameters
    al = 0.5 * (gamma + 1.0)
    ar = 0.5 * (gamma + 1.0)

    csl, csr, umin, umax, plmin, prmin, bl, br = declare('double', 8)
    a, b, c, d, pstar, ustar, dd = declare('double', 7)
    # Lagrangian sound speeds
    csl = sqrt(gamma * pl * rhol)
    csr = sqrt(gamma * pr * rhor)

    umin = ur - 0.5 * csr / ar
    umax = ul + 0.5 * csl / al

    plmin = pl - 0.25 * rhol * csl * csl / al
    prmin = pr - 0.25 * rhor * csr * csr / ar

    bl = rhol * al
    br = rhor * ar

    a = (br - bl) * (prmin - plmin)
    b = br * umin * umin - bl * umax * umax
    c = br * umin - bl * umax
    d = br * bl * (umin - umax) * (umin - umax)

    # Case A : ustar - umin > 0
    dd = sqrt(max(0.0, d - a))
    ustar = (b + prmin - plmin) / (c - SIGN(dd, umax - umin))

    if (((ustar - umin) >= 0.0) and ((ustar - umax) <= 0)):
        pstar = 0.5 * (plmin + prmin + br * abs(ustar - umin) *
                       (ustar - umin) - bl * abs(ustar - umax) *
                       (ustar - umax))
        pstar = max(pstar, 0.0)
        result[0] = pstar
        result[1] = ustar
        return 0

    # Case B: ustar - umin < 0, ustar - umax > 0
    dd = sqrt(max(0.0, d + a))
    ustar = (b - prmin + plmin) / (c - SIGN(dd, umax - umin))
    if (((ustar - umin) <= 0.0) and ((ustar - umax) >= 0.0)):
        pstar = 0.5 * (plmin + prmin + br * abs(ustar - umin) *
                       (ustar - umin) - bl * abs(ustar - umax) *
                       (ustar - umax))
        pstar = max(pstar, 0.0)
        result[0] = pstar
        result[1] = ustar
        return 0

    a = (bl + br) * (plmin - prmin)
    b = bl * umax + br * umin
    c = 1. / (bl + br)

    # Case C : ustar-umin >0, ustar-umax > 0
    dd = sqrt(max(0.0, a - d))
    ustar = (b + dd) * c
    if (((ustar - umin) >= 0) and ((ustar - umax) >= 0.0)):
        pstar = 0.5 * (plmin + prmin + br * abs(ustar - umin) *
                       (ustar - umin) - bl * abs(ustar - umax) *
                       (ustar - umax))
        pstar = max(pstar, 0.0)
        result[0] = pstar
        result[1] = ustar
        return 0

    # Case D: ustar - umin < 0, ustar - umax < 0
    dd = sqrt(max(0.0, -a - d))
    ustar = (b - dd) * c
    pstar = 0.5 * (plmin + prmin + br * abs(ustar - umin) *
                   (ustar - umin) - bl * abs(ustar - umax) * (ustar - umax))
    pstar = max(pstar, 0.0)
    result[0] = pstar
    result[1] = ustar
    return 0
Пример #19
0
 def reduce(self, dst):
     indices = declare('object')
     indices = numpy.where(dst.merge > 0)[0]
     if numpy.any(indices):
         dst.remove_particles(indices)
Пример #20
0
def van_leer(rhol=0.0,
             rhor=1.0,
             pl=0.0,
             pr=1.0,
             ul=0.0,
             ur=1.0,
             gamma=1.4,
             niter=20,
             tol=1e-6,
             result=[0.0, 0.0]):
    """Van Leer Riemann solver.


    Parameters
    ----------
    rhol, rhor: double: left and right density.
    pl, pr: double: left and right pressure.
    ul, ur: double: left and right speed.
    gamma: double: Ratio of specific heats.
    niter: int: Max number of iterations to try for iterative schemes.
    tol: double: Error tolerance for convergence.

    result: list: List of length 2. Result will contain computed pstar, ustar

    Returns
    -------

    Returns 0 if the value is computed correctly else 1 if the iterations
    either did not converge or if there is an error.

    """
    if ((rhol < 0) or (rhor < 0) or (pl < 0) or (pr < 0)):
        result[0] = 0.0
        result[1] = 0.0
        return 1

    # local variables
    gamma1, gamma2 = declare('double', 2)
    Vl, Vr, cl, cr, zl, zr, wl, wr = declare('double', 8)
    ustar_l, ustar_r, pstar, pstar_old = declare('double', 4)
    converged, iteration = declare('int', 2)

    gamma2 = 1.0 + gamma
    gamma1 = 0.5 * gamma2 / gamma
    smallp = 1e-25

    # initialize variables
    wl = 0.0
    wr = 0.0
    zl = 0.0
    zr = 0.0

    # specific volumes
    Vl = 1. / rhol
    Vr = 1. / rhor

    # Lagrangean sound speeds
    cl = sqrt(gamma * pl * rhol)
    cr = sqrt(gamma * pr * rhor)

    # Initial guess for pstar
    pstar = pr - pl - cr * (ur - ul)
    pstar = pl + pstar * cl / (cl + cr)

    pstar = max(pstar, smallp)

    # Now Iterate using NR to obtain the star values
    iteration = 0
    while iteration < niter:
        pstar_old = pstar

        wl = 1.0 + gamma1 * (pstar - pl) / pl
        wl = cl * sqrt(wl)

        wr = 1.0 + gamma1 * (pstar - pr) / pr
        wr = cr * sqrt(wr)

        zl = 4.0 * Vl * wl * wl
        zl = -zl * wl / (zl - gamma2 * (pstar - pl))

        zr = 4.0 * Vr * wr * wr
        zr = zr * wr / (zr - gamma2 * (pstar - pr))

        ustar_l = ul - (pstar - pl) / wl
        ustar_r = ur + (pstar - pr) / wr

        pstar = pstar + (ustar_r - ustar_l) * (zl * zr) / (zr - zl)
        pstar = max(smallp, pstar)

        converged = (abs(pstar - pstar_old) / pstar < tol)
        if converged:
            break

        iteration += 1

    # calculate the averaged ustar
    ustar_l = ul - (pstar - pl) / wl
    ustar_r = ur + (pstar - pr) / wr
    ustar = 0.5 * (ustar_l + ustar_r)

    result[0] = pstar
    result[1] = ustar
    if converged:
        return 0
    else:
        return 1
Пример #21
0
 def func(self, d_idx, d_x=[0.0, 0.0]):
     mat = declare('matrix((2,2))')
     mat[0][0] = d_x[d_idx]
     vec, vec1 = declare('matrix(3)', 2)
     vec[0] = d_x[d_idx]
Пример #22
0
def hllc_ball(rhol=0.0,
              rhor=1.0,
              pl=0.0,
              pr=1.0,
              ul=0.0,
              ur=1.0,
              gamma=1.4,
              niter=20,
              tol=1e-6,
              result=[0.0, 0.0]):
    """Harten-Lax-van Leer-Contact approximate Riemann solver for the Euler
    equations to determine the intermediate states.

    Parameters
    ----------
    rhol, rhor: double: left and right density.
    pl, pr: double: left and right pressure.
    ul, ur: double: left and right speed.
    gamma: double: Ratio of specific heats.
    niter: int: Max number of iterations to try for iterative schemes.
    tol: double: Error tolerance for convergence.

    result: list: List of length 2. Result will contain computed pstar, ustar

    Returns
    -------

    Returns 0 if the value is computed correctly else 1 if the iterations
    either did not converge or if there is an error.

    """
    gamma1, csl, csr, cslr, rholr, pstar, ustar = declare('double', 7)
    Sl, Sr, ql, qr, pstar_l, pstar_r = declare('double', 6)

    gamma1 = 0.5 * (gamma + 1.0) / gamma

    # sound speeds in the undisturbed fluid
    csl = sqrt(gamma * pl / rhol)
    csr = sqrt(gamma * pr / rhor)
    cslr = 0.5 * (csl + csr)

    # averaged density
    rholr = 0.5 * (rhol + rhor)

    # provisional intermediate pressure
    pstar = 0.5 * (pl + pr - rholr * cslr * (ur - ul))

    # Wave speed estimates (ustar is taken as the intermediate velocity)
    ustar = 0.5 * (ul + ur - 1. / (rholr * cslr) * (pr - pl))

    Hl = pstar / pl
    Hr = pstar / pr

    ql = 1.0
    if (Hl > 1):
        ql = sqrt(1 + gamma1 * (Hl - 1.0))

    qr = 1.0
    if (Hr > 1):
        qr = sqrt(1 + gamma1 * (Hr - 1.0))

    Sl = ul - csl * ql
    Sr = ur + csr * qr

    # compute the intermediate pressure
    pstar_l = pl + rhol * (ul - Sl) * (ul - ustar)
    pstar_r = pr + rhor * (ur - Sr) * (ur - ustar)

    pstar = 0.5 * (pstar_l + pstar_r)

    # return intermediate states
    result[0] = pstar
    result[1] = ustar
    return 0
Пример #23
0
    def test_declare(self):
        self.assertEqual(declare('int'), 0)
        self.assertEqual(declare('long'), 0)
        self.assertEqual(declare('double'), 0.0)
        self.assertEqual(declare('float'), 0.0)

        self.assertEqual(declare('int', 2), (0, 0))
        self.assertEqual(declare('long', 3), (0, 0, 0))
        self.assertEqual(declare('double', 2), (0.0, 0.0))
        self.assertEqual(declare('float', 3), (0.0, 0.0, 0.0))

        res = declare('matrix(3)')
        self.assertTrue(np.all(res == np.zeros(3)))
        res = declare('matrix(3)', 3)
        for i in range(3):
            self.assertTrue(np.all(res[0] == np.zeros(3)))
        res = declare('matrix((3,))')
        self.assertTrue(np.all(res == np.zeros(3)))
        res = declare('matrix((3, 3))')
        self.assertTrue(np.all(res == np.zeros((3, 3))))