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]
def velocity(i, x, y, gamma, u, v, nv): j = declare('int') tmp = declare('matrix(2)') xi = x[i] yi = y[i] u[i] = 0.0 v[i] = 0.0 for j in range(nv): point_vortex(xi, yi, x[j], y[j], gamma[j], tmp) u[i] += tmp[0] v[i] += tmp[1]
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
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 = abs(xij[0]) sum += s_m[s_idx]*KERNEL.kernel(xij, rij, s_h[s_idx]) d_rho[d_idx] += sum
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, SPH_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 = SPH_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 = SPH_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
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
def align_output_expr(i, item, prev_item, last_item, tag_arr, new_indices, num_particles, stride): t, idx, j_s = declare('int', 3) t = last_item + i - prev_item idx = t if tag_arr[i] else prev_item for j_s in range(stride): new_indices[stride * idx + j_s] = stride * i + j_s
def monotonicity_min(_x1=0.0, _x2=0.0, _x3=0.0): x1, x2, x3, _min = declare('double', 4) x1 = 2.0 * abs(_x1) x2 = abs(_x2) x3 = 2.0 * abs(_x3) sx1 = sgn(_x1) sx2 = sgn(_x2) sx3 = sgn(_x3) if ((sx1 != sx2) or (sx2 != sx3)): return 0.0 else: if x2 < x1: if x3 < x2: _min = x3 else: _min = x2 else: if x3 < x1: _min = x3 else: _min = x1 return sx1 * _min
def align_output_expr(i, item, prev_item, last_item, tag_arr, new_indices, num_particles, num_real_particles): t, idx = declare('int', 2) t = last_item + i - prev_item idx = t if tag_arr[i] else prev_item new_indices[idx] = i if i == num_particles - 1: num_real_particles[0] = last_item
def reduce(self, dst, t, dt): 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
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, SPH_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 = SPH_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
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
def stage1(self, d_idx, d_x, d_y, d_z, d_u, d_v, d_w, d_fx, d_fy, d_fz, d_torx, d_tory, d_torz, d_wx, d_wy, d_wz, d_m_inverse, d_I_inverse, d_total_dem_entities, d_total_tng_contacts, d_limit, d_tng_x, d_tng_y, d_tng_z, d_tng_x0, d_tng_y0, d_tng_z0, d_vtx, d_vty, d_vtz, dt): dtb2 = dt / 2. i = declare('int') j = declare('int') p = declare('int') q = declare('int') idx_total_ctcs_with_entity_i = declare('int') # loop over all entites for i in range(0, d_total_dem_entities[0]): idx_total_ctcs_with_entity_i = ( d_total_tng_contacts[d_total_dem_entities[0] * d_idx + i]) # particle idx contacts with entity i has range of indices # and the first index would be p = d_idx * d_total_dem_entities[0] * d_limit[0] + i * d_limit[0] q = p + idx_total_ctcs_with_entity_i # loop over all the contacts of particle d_idx with entity i for j in range(p, q): d_tng_x[j] = d_tng_x0[j] + d_vtx[j] * dtb2 d_tng_y[j] = d_tng_y0[j] + d_vty[j] * dtb2 d_tng_z[j] = d_tng_z0[j] + d_vtz[j] * dtb2 d_x[d_idx] = d_x[d_idx] + dtb2 * d_u[d_idx] d_y[d_idx] = d_y[d_idx] + dtb2 * d_v[d_idx] d_z[d_idx] = d_z[d_idx] + dtb2 * d_w[d_idx] d_u[d_idx] = d_u[d_idx] + dtb2 * d_fx[d_idx] * d_m_inverse[d_idx] d_v[d_idx] = d_v[d_idx] + dtb2 * d_fy[d_idx] * d_m_inverse[d_idx] d_w[d_idx] = d_w[d_idx] + dtb2 * d_fz[d_idx] * d_m_inverse[d_idx] d_wx[d_idx] = d_wx[d_idx] + dtb2 * d_torx[d_idx] * d_I_inverse[d_idx] d_wy[d_idx] = d_wy[d_idx] + dtb2 * d_tory[d_idx] * d_I_inverse[d_idx] d_wz[d_idx] = d_wz[d_idx] + dtb2 * d_torz[d_idx] * d_I_inverse[d_idx]
def velocity(x, y, gamma, u, v, xc, yc, gc, nv): i, gid, nb = declare('int', 3) j, ti, nt, jb = declare('int', 4) ti = LID_0 nt = LDIM_0 gid = GID_0 i = gid*nt + ti idx = declare('int') tmp = declare('matrix(2)') uj, vj = declare('double', 2) nb = GDIM_0 if i < nv: xi = x[i] yi = y[i] uj = 0.0 vj = 0.0 for jb in range(nb): idx = jb*nt + ti if idx < nv: xc[ti] = x[idx] yc[ti] = y[idx] gc[ti] = gamma[idx] else: gc[ti] = 0.0 local_barrier() if i < nv: for j in range(nt): point_vortex(xi, yi, xc[j], yc[j], gc[j], tmp) uj += tmp[0] vj += tmp[1] local_barrier() if i < nv: u[i] = uj v[i] = vj
def interpolate(self, hi=0.0, hj=0.0, rhoi=0.0, rhoj=0.0, sij=0.0, gri_eij=0.0, grj_eij=0.0, result=[0.0, 0.0, 0.0]): """Interpolation for the specific volume integrals in GSPH. Parameters: ----------- hi, hj : double Particle smoothing length (scale?) at i (right) and j (left) rhoi, rhoj : double Particle densities at i (right) and j (left) sij : double Particle separation in the local coordinate system (si - sj) gri_eij, grj_eij : double Gradient of density at the particles in the global coordinate system dot product with eij Notes: ------ The interpolation scheme determines the form of the 'specific volume' contributions Vij^2 in the GSPH equations. The simplest Delta or point interpolation uses Vi^2 = 1./rho_i^2. The more involved linear or cubic spline interpolations are defined in the GSPH references. Most recent papers on GSPH typically assume the interface to be located midway between the particles. In the local coordinate system, this corresponds to sij = 0. From I02, it seems the definition is only really valid for the constant smoothing length case. Set interface_zero to False if you want to include this term """ Vi = 1. / rhoi Vj = 1. / rhoj Vip = -1. / (rhoi * rhoi) * gri_eij Vjp = -1. / (rhoj * rhoj) * grj_eij aij, bij, cij, dij, hij, vij, vij_i, vij_j = declare('double', 8) hij = 0.5 * (hi + hj) sstar = self.sstar # simplest delta or point interpolation if self.interpolation == 0: vij_i = 1. / (rhoi * rhoi) vij_j = 1. / (rhoj * rhoj) # linear interpolation elif self.interpolation == 1: # avoid singularities if sij < 1e-8: cij = 0.0 else: cij = (Vi - Vj) / sij dij = 0.5 * (Vi + Vj) vij_i = 0.25 * hi * hi * cij * cij + dij * dij vij_j = 0.25 * hj * hj * cij * cij + dij * dij # approximate value for the interface location when using # variable smoothing lengths if not self.interface_zero: vij = 0.5 * (vij_i + vij_j) sstar = 0.5 * hij * hij * cij * dij / vij # cubic spline interpolation elif self.interpolation == 2: if sij < 1e-8: aij = bij = cij = 0.0 dij = 0.5 * Vi + Vj else: aij = -2.0 * (Vi - Vj) / (sij * sij * sij) + (Vip + Vjp) / (sij * sij) bij = 0.5 * (Vip - Vjp) / sij cij = 1.5 * (Vi - Vj) / sij - 0.25 * (Vip + Vjp) dij = 0.5 * (Vi + Vj) - 0.125 * (Vip - Vjp) * sij vij_i = ((15.0) / (64.0) * hi**6 * aij * aij + (3.0) / (16.0) * hi**4 * (2 * aij * cij + bij * bij) + 0.25 * hi * hi * (2 * bij * dij + cij * cij) + dij * dij) vij_j = ((15.0) / (64.0) * hj**6 * aij * aij + (3.0) / (16.0) * hj**4 * (2 * aij * cij + bij * bij) + 0.25 * hj * hj * (2 * bij * dij + cij * cij) + dij * dij) else: printf("Unknown interpolation type") result[0] = vij_i result[1] = vij_j result[2] = sstar
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
def loop(self, d_idx, d_m, d_wx, d_wy, d_wz, d_fx, d_fy, d_fz, d_torx, d_tory, d_torz, d_tng_x, d_tng_y, d_tng_z, d_tng_x0, d_tng_y0, d_tng_z0, d_tng_idx, d_total_tng_contacts, d_dem_id, d_limit, d_vtx, d_vty, d_vtz, d_total_dem_entities, VIJ, XIJ, RIJ, d_R, s_idx, s_R, s_wx, s_wy, s_wz, s_dem_id, dt): overlap = -1. print(overlap, "is good") # check if we are not dealing with the same particle if RIJ > 0: overlap = d_R[d_idx] + s_R[s_idx] - RIJ # d_idx has a range of tracking indices with source dem_id. # starting index is p p = declare('int') # ending index is q -1 q = declare('int') # total number of contacts of particle i in destination # with source entity tot_ctcs = declare('int') tot_ctcs = (d_total_tng_contacts[d_idx * d_total_dem_entities[0] + s_dem_id[0]]) p = (d_idx * d_total_dem_entities[0] * d_limit[0] + s_dem_id[0] * d_limit[0]) q = p + tot_ctcs i = declare('int') found_at = declare('int') found = declare('int') # check if the particle is in the tracking list # if so, then save the location at found_at found = 0 for i in range(p, q): if s_idx == d_tng_idx[i]: found_at = i found = 1 break # ---------- force computation starts ------------ # if particles are not overlapping if overlap <= 0: if found == 1: # make its tangential displacement to be zero d_tng_x0[found_at] = 0. d_tng_y0[found_at] = 0. d_tng_z0[found_at] = 0. d_tng_x[found_at] = 0. d_tng_y[found_at] = 0. d_tng_z[found_at] = 0. d_vtx[found_at] = 0. d_vty[found_at] = 0. d_vtz[found_at] = 0. # if particles are in contact else: # normal vector passing from s_idx to d_idx, i.e., j to i rinv = 1.0 / RIJ nx = XIJ[0] * rinv ny = XIJ[1] * rinv nz = XIJ[2] * rinv # ---- Realtive velocity computation (Eq 10)---- # relative velocity of particle d_idx w.r.t particle s_idx at # contact point. The velocity different provided by PySPH is # only between translational velocities, but we need to consider # rotational velocities also. # Distance till contact point a_i = d_R[d_idx] - overlap / 2. a_j = s_R[s_idx] - overlap / 2. vr_x = VIJ[0] + (a_i * (ny * d_wz[d_idx] - nz * d_wy[d_idx]) + a_j * (ny * s_wz[s_idx] - nz * s_wy[s_idx])) vr_y = VIJ[1] + (a_i * (nz * d_wx[d_idx] - nx * d_wz[d_idx]) + a_j * (nz * s_wx[s_idx] - nx * s_wz[s_idx])) vr_z = VIJ[2] + (a_i * (nx * d_wy[d_idx] - ny * d_wx[d_idx]) + a_j * (nx * s_wy[s_idx] - ny * s_wx[s_idx])) # normal velocity magnitude vr_dot_nij = vr_x * nx + vr_y * ny + vr_z * nz vn_x = vr_dot_nij * nx vn_y = vr_dot_nij * ny vn_z = vr_dot_nij * nz # tngential velocity vt_x = vr_x - vn_x vt_y = vr_y - vn_y vt_z = vr_z - vn_z # tangential velocity magnitude # vt_magn = (vt_x * vt_x + vt_y * vt_y + vt_z * vt_z)**(1. / 2.) # normal force fn_x = self.kn * overlap * nx - self.eta_n * vn_x fn_y = self.kn * overlap * ny - self.eta_n * vn_y fn_z = self.kn * overlap * nz - self.eta_n * vn_z # ------------- tangential force computation ---------------- # if the particle is not been tracked then assign an index in # tracking history. if found == 0: found_at = q d_tng_idx[found_at] = s_idx d_total_tng_contacts[d_idx * d_total_dem_entities[0] + s_dem_id[0]] += 1 # compute and set the tangential acceleration for the current # time step d_vtx[found_at] = vt_x d_vty[found_at] = vt_y d_vtz[found_at] = vt_z else: # rotate the spring for current plane tng_dot_nij = (d_tng_x[found_at] * nx + d_tng_y[found_at] * ny + d_tng_z[found_at] * nz) d_tng_x[found_at] -= tng_dot_nij * nx d_tng_y[found_at] -= tng_dot_nij * ny d_tng_z[found_at] -= tng_dot_nij * nz # compute and set the tangential acceleration for the current # time step d_vtx[found_at] = vt_x d_vty[found_at] = vt_y d_vtz[found_at] = vt_z # find the tangential force from the tangential displacement # and tangential velocity (eq 18 Luding) ft0_x = -self.kt * d_tng_x[found_at] - self.eta_t * vt_x ft0_y = -self.kt * d_tng_y[found_at] - self.eta_t * vt_y ft0_z = -self.kt * d_tng_z[found_at] - self.eta_t * vt_z # (*) check against Coulomb criterion # Tangential force magnitude due to displacement ft0_magn = (ft0_x * ft0_x + ft0_y * ft0_y + ft0_z * ft0_z)**(0.5) fn_magn = (fn_x * fn_x + fn_y * fn_y + fn_z * fn_z)**(0.5) # we have to compare with static friction, so # this mu has to be static friction coefficient fn_mu = self.mu * fn_magn # if the tangential force magnitude is zero, then do nothing, # else do following if ft0_magn != 0.: # compare tangential force with the static friction if ft0_magn >= fn_mu: # rescale the tangential displacement tx = ft0_x / ft0_magn ty = ft0_y / ft0_magn tz = ft0_z / ft0_magn d_tng_x[found_at] = -self.kt_1 * ( fn_mu * tx + self.eta_t * vt_x) d_tng_y[found_at] = -self.kt_1 * ( fn_mu * ty + self.eta_t * vt_y) d_tng_z[found_at] = -self.kt_1 * ( fn_mu * tz + self.eta_t * vt_z) # set the tangential force to static friction # from Coulomb criterion ft0_x = fn_mu * tx ft0_y = fn_mu * ty ft0_z = fn_mu * tz d_fx[d_idx] += fn_x + ft0_x d_fy[d_idx] += fn_y + ft0_y d_fz[d_idx] += fn_z + ft0_z d_torx[d_idx] += (ft0_y * nz - ft0_z * ny) * a_i d_tory[d_idx] += (ft0_z * nx - ft0_x * nz) * a_i d_torz[d_idx] += (ft0_x * ny - ft0_y * nx) * a_i
def initialize_pair(self, d_idx, d_x, d_y, d_z, d_R, d_total_dem_entities, d_total_tng_contacts, d_tng_idx, d_limit, d_tng_x, d_tng_y, d_tng_z, s_x, s_y, s_z, s_R): i = declare('int') p = declare('int') count = declare('int') k = declare('int') xij = declare('matrix(3)') last_idx_tmp = declare('int') sidx = declare('int') rij = 0.0 idx_total_ctcs_with_entity_i = declare('int') # loop over all entites for i in range(0, d_total_dem_entities[0]): idx_total_ctcs_with_entity_i = ( d_total_tng_contacts[d_total_dem_entities[0] * d_idx + i]) # particle idx contacts with entity i has range of indices # and the first index would be p = d_idx * d_total_dem_entities[0] * d_limit[0] + i * d_limit[0] last_idx_tmp = p + idx_total_ctcs_with_entity_i - 1 k = p count = 0 # loop over all the contacts of particle d_idx with entity i while count < idx_total_ctcs_with_entity_i: # The index of the particle with which # d_idx in contact is sidx = d_tng_idx[k] if sidx == -1: break else: xij[0] = d_x[d_idx] - s_x[sidx] xij[1] = d_y[d_idx] - s_y[sidx] xij[2] = d_z[d_idx] - s_z[sidx] rij = sqrt(xij[0] * xij[0] + xij[1] * xij[1] + xij[2] * xij[2]) overlap = d_R[d_idx] + s_R[sidx] - rij if overlap < 0.: # if the swap index is the current index then # simply make it to null contact. if k == last_idx_tmp: d_tng_idx[k] = -1 d_tng_x[k] = 0. d_tng_y[k] = 0. d_tng_z[k] = 0. else: # swap the current tracking index with the final # contact index d_tng_idx[k] = d_tng_idx[last_idx_tmp] d_tng_idx[last_idx_tmp] = -1 # swap tangential x displacement d_tng_x[k] = d_tng_x[last_idx_tmp] d_tng_x[last_idx_tmp] = 0. # swap tangential y displacement d_tng_y[k] = d_tng_x[last_idx_tmp] d_tng_y[last_idx_tmp] = 0. # swap tangential z displacement d_tng_z[k] = d_tng_z[last_idx_tmp] d_tng_z[last_idx_tmp] = 0. # decrease the last_idx_tmp, since we swapped it to # -1 last_idx_tmp -= 1 else: k = k + 1 count += 1
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
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
def fill_indices_elwise(i, indices, new_indices, size, stride): tmp_idx, j_s = declare('unsigned int', 2) for j_s in range(stride): tmp_idx = i * stride + j_s if tmp_idx < size: new_indices[tmp_idx] = indices[i] * stride + j_s
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
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
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
def loop(self, d_idx, d_m, d_h, d_rho, d_cs, d_div, d_p, d_e, d_grhox, d_grhoy, d_grhoz, d_u, d_v, d_w, d_px, d_py, d_pz, d_ux, d_uy, d_uz, d_vx, d_vy, d_vz, d_wx, d_wy, d_wz, d_au, d_av, d_aw, d_ae, s_idx, s_rho, s_m, s_h, s_cs, s_div, s_p, s_e, s_grhox, s_grhoy, s_grhoz, s_u, s_v, s_w, s_px, s_py, s_pz, s_ux, s_uy, s_uz, s_vx, s_vy, s_vz, s_wx, s_wy, s_wz, XIJ, DWIJ, DWI, DWJ, RIJ, RHOIJ, EPS, dt, t): blending_factor = exp(-self.blend_alpha * t / self.tf) g1 = self.g1 g2 = self.g2 hi = d_h[d_idx] hj = s_h[s_idx] eij = declare('matrix(3)') if RIJ < 1e-14: eij[0] = 0.0 eij[1] = 0.0 eij[2] = 0.0 sij = 1.0 / (RIJ + EPS) else: eij[0] = XIJ[0] / RIJ eij[1] = XIJ[1] / RIJ eij[2] = XIJ[2] / RIJ sij = 1.0 / RIJ vl = s_u[s_idx] * eij[0] + s_v[s_idx] * eij[1] + s_w[s_idx] * eij[2] vr = d_u[d_idx] * eij[0] + d_v[d_idx] * eij[1] + d_w[d_idx] * eij[2] Hi = g1 * hi * d_cs[d_idx] + g2 * hi * hi * (abs(d_div[d_idx]) - d_div[d_idx]) grhoi_dot_eij = (d_grhox[d_idx] * eij[0] + d_grhoy[d_idx] * eij[1] + d_grhoz[d_idx] * eij[2]) grhoj_dot_eij = (s_grhox[s_idx] * eij[0] + s_grhoy[s_idx] * eij[1] + s_grhoz[s_idx] * eij[2]) sp_vol = declare('matrix(3)') self.interpolate(hi, hj, d_rho[d_idx], s_rho[s_idx], RIJ, grhoi_dot_eij, grhoj_dot_eij, sp_vol) vij_i = sp_vol[0] vij_j = sp_vol[1] sstar = sp_vol[2] # Gradients in the local coordinate system rsi = (d_grhox[d_idx] * eij[0] + d_grhoy[d_idx] * eij[1] + d_grhoz[d_idx] * eij[2]) psi = (d_px[d_idx] * eij[0] + d_py[d_idx] * eij[1] + d_pz[d_idx] * eij[2]) vsi = (eij[0] * eij[0] * d_ux[d_idx] + eij[0] * eij[1] * (d_uy[d_idx] + d_vx[d_idx]) + eij[0] * eij[2] * (d_uz[d_idx] + d_wx[d_idx]) + eij[1] * eij[1] * d_vy[d_idx] + eij[1] * eij[2] * (d_vz[d_idx] + d_wy[d_idx]) + eij[2] * eij[2] * d_wy[d_idx]) rsj = (s_grhox[s_idx] * eij[0] + s_grhoy[s_idx] * eij[1] + s_grhoz[s_idx] * eij[2]) psj = (s_px[s_idx] * eij[0] + s_py[s_idx] * eij[1] + s_pz[s_idx] * eij[2]) vsj = (eij[0] * eij[0] * s_ux[s_idx] + eij[0] * eij[1] * (s_uy[s_idx] + s_vx[s_idx]) + eij[0] * eij[2] * (s_uz[s_idx] + s_wx[s_idx]) + eij[1] * eij[1] * s_vy[s_idx] + eij[1] * eij[2] * (s_vz[s_idx] + s_wy[s_idx]) + eij[2] * eij[2] * s_wy[s_idx]) csi = d_cs[d_idx] csj = s_cs[s_idx] rhoi = d_rho[d_idx] rhoj = s_rho[s_idx] pi = d_p[d_idx] pj = s_p[s_idx] if self.monotonicity == 0: # First order scheme rsi = 0.0 rsj = 0.0 psi = 0.0 psj = 0.0 vsi = 0.0 vsj = 0.0 if self.monotonicity == 1: # I02 algorithm if (vsi * vsj) < 0: vsi = 0. vsj = 0. # default to first order near a shock if (min(csi, csj) < 3.0 * (vl - vr)): rsi = 0. rsj = 0. psi = 0. psj = 0. vsi = 0. vsj = 0. if self.monotonicity == 2 and RIJ > 1e-14: # IwIn algorithm qijr = rhoi - rhoj qijp = pi - pj qiju = vr - vl delr = rsi * RIJ delrp = 2 * delr - qijr delp = psi * RIJ delpp = 2 * delp - qijp delv = vsi * RIJ delvp = 2 * delv - qiju # corrected values for i rsi = monotonicity_min(qijr, delr, delrp) / RIJ psi = monotonicity_min(qijp, delp, delpp) / RIJ vsi = monotonicity_min(qiju, delv, delvp) / RIJ delr = rsj * RIJ delrp = 2 * delr - qijr delp = psj * RIJ delpp = 2 * delp - qijp delv = vsj * RIJ delvp = 2 * delv - qiju # corrected values for j rsj = monotonicity_min(qijr, delr, delrp) / RIJ psj = monotonicity_min(qijp, delp, delpp) / RIJ vsj = monotonicity_min(qiju, delv, delvp) / RIJ elif self.monotonicity == 2 and RIJ < 1e-14: # IwIn algorithm rsi = 0.0 rsj = 0.0 psi = 0.0 psj = 0.0 vsi = 0.0 vsj = 0.0 # Input to the riemann solver sstar *= 2.0 # left and right density rhol = rhoj + 0.5 * rsj * RIJ * (1.0 - csj * dt * sij + sstar) rhor = rhoi - 0.5 * rsi * RIJ * (1.0 - csi * dt * sij + sstar) # corrected density if rhol < 0: rhol = rhoj if rhor < 0: rhor = rhoi # left and right pressure pl = pj + 0.5 * psj * RIJ * (1.0 - csj * dt * sij + sstar) pr = pi - 0.5 * psi * RIJ * (1.0 - csi * dt * sij + sstar) # corrected pressure if pl < 0: pl = pj if pr < 0: pr = pi # left and right velocity ul = vl + 0.5 * vsj * RIJ * (1.0 - csj * dt * sij + sstar) ur = vr - 0.5 * vsi * RIJ * (1.0 - csi * dt * sij + sstar) # Intermediate state from the Riemann solver result = declare('matrix(2)') riemann_solve(self.rsolver, rhol, rhor, pl, pr, ul, ur, self.gamma, self.niter, self.tol, result) pstar = result[0] ustar = result[1] # blend of two intermediate states if self.hybrid: riemann_solve(10, rhoj, rhoi, pl, pr, vl, vr, self.gamma, self.niter, self.tol, result) pstar2 = result[0] ustar2 = result[1] ustar = ustar + blending_factor * (ustar2 - ustar) pstar = pstar + blending_factor * (pstar2 - pstar) # three dimensional velocity (70) vstar = declare('matrix(3)') vstar[0] = ustar * eij[0] vstar[1] = ustar * eij[1] vstar[2] = ustar * eij[2] # velocity accelerations mj = s_m[s_idx] d_au[d_idx] += -mj * pstar * (vij_i * DWI[0] + vij_j * DWJ[0]) d_av[d_idx] += -mj * pstar * (vij_i * DWI[1] + vij_j * DWJ[1]) d_aw[d_idx] += -mj * pstar * (vij_i * DWI[2] + vij_j * DWJ[2]) # contribution to the thermal energy term # (85). The contribution due to \dot{x}^*_i will # be added in the integrator. vstardotdwi = vstar[0] * DWI[0] + vstar[1] * DWI[1] + vstar[2] * DWI[2] vstardotdwj = vstar[0] * DWJ[0] + vstar[1] * DWJ[1] + vstar[2] * DWJ[2] d_ae[d_idx] += -mj * pstar * (vij_i * vstardotdwi + vij_j * vstardotdwj) # artificial thermal conduction terms if self.thermal_conduction: divj = s_div[s_idx] Hj = g1 * hj * csj + g2 * hj * hj * (abs(divj) - divj) Hij = (Hi + Hj) * (d_e[d_idx] - s_e[s_idx]) Hij /= (RHOIJ * (RIJ * RIJ + EPS)) d_ae[d_idx] += mj * Hij * (XIJ[0] * DWIJ[0] + XIJ[1] * DWIJ[1] + XIJ[2] * DWIJ[2])