def doit(nx, ny): # test the multigrid solver # create the multigrid object a = MG.CellCenterMG2d(nx, ny, xl_BC_type="dirichlet", yl_BC_type="dirichlet", xr_BC_type="dirichlet", yr_BC_type="dirichlet", verbose=0, nsmooth=5, nsmooth_bottom=10, vis=1, true_function=true, vis_title=r"$u_{xx} + u_{yy} = -2[(1-6x^2)y^2(1-y^2) + (1-6y^2)x^2(1-x^2)]$") plt.ion() plt.figure(num=1, figsize=(12.8, 7.2), dpi=100, facecolor='w') # initialize the solution to 0 init = a.soln_grid.scratch_array() a.init_solution(init) # initialize the RHS using the function f rhs = f(a.x2d, a.y2d) a.init_RHS(rhs) # solve to a relative tolerance of 1.e-11 a.solve(rtol=1.e-11) # alternately, we can just use smoothing by uncommenting the following # a.smooth(a.nlevels-1,50000) # get the solution v = a.get_solution() # compute the error from the analytic solution b = true(a.x2d, a.y2d) e = v - b print(" L2 error from true solution = %g\n rel. err from previous cycle = %g\n num. cycles = %d" % (a.soln_grid.norm(e), a.relative_error, a.num_cycles)) # plot it # plt.figure(num=1, figsize=(2.10,2.10), dpi=100, facecolor='w') plt.figure(num=1, figsize=(5.0, 5.0), dpi=100, facecolor='w') plt.imshow(np.transpose(v[a.ilo:a.ihi+1, a.jlo:a.jhi+1]), interpolation="nearest", origin="lower", extent=[a.xmin, a.xmax, a.ymin, a.ymax]) # plt.axis("off") # plt.subplots_adjust(bottom=0.0, top=1.0, left=0.0, right=1.0) plt.xlabel("x") plt.ylabel("y") plt.savefig("mg_test.png") # store the output for later comparison my_data = a.get_solution_object() my_data.write("mg_test")
def evolve(self): """ Diffusion through dt using C-N implicit solve with multigrid """ self.cc_data.fill_BC_all() phi = self.cc_data.get_var("phi") myg = self.cc_data.grid # diffusion coefficient k = self.rp.get_param("diffusion.k") # setup the MG object -- we want to solve a Helmholtz equation # equation of the form: # (alpha - beta L) phi = f # # with alpha = 1 # beta = (dt/2) k # f = phi + (dt/2) k L phi # # this is the form that arises with a Crank-Nicolson discretization # of the diffusion equation. mg = MG.CellCenterMG2d(myg.nx, myg.ny, xmin=myg.xmin, xmax=myg.xmax, ymin=myg.ymin, ymax=myg.ymax, xl_BC_type=self.cc_data.BCs['phi'].xlb, xr_BC_type=self.cc_data.BCs['phi'].xrb, yl_BC_type=self.cc_data.BCs['phi'].ylb, yr_BC_type=self.cc_data.BCs['phi'].yrb, alpha=1.0, beta=0.5 * self.dt * k, verbose=0) # form the RHS: f = phi + (dt/2) k L phi (where L is the Laplacian) f = mg.soln_grid.scratch_array() f.v()[:, :] = phi.v() + 0.5 * self.dt * k * ( (phi.ip(1) + phi.ip(-1) - 2.0 * phi.v()) / myg.dx**2 + (phi.jp(1) + phi.jp(-1) - 2.0 * phi.v()) / myg.dy**2) mg.init_RHS(f) # initial guess is zeros mg.init_zeros() # solve the MG problem for the updated phi mg.solve(rtol=1.e-10) #mg.smooth(mg.nlevels-1,100) # update the solution phi.v()[:, :] = mg.get_solution().v() # increment the time self.cc_data.t += self.dt self.n += 1
def test_mg_gradient(): a = MG.CellCenterMG2d(8, 8, ng=1, xmax=8, ymax=8) s = a.soln_grid.scratch_array() s.v()[:,:] = np.fromfunction(lambda i, j: i*(s.g.nx-i-1)*j*(s.g.ny-j-1), (s.g.nx, s.g.ny)) a.init_solution(s) a.grids[a.nlevels-1].fill_BC("v") gx, gy = a.get_solution_gradient() assert_array_equal(gx[:,gx.g.jc], np.array([0., 36., 60., 36., 12., -12., -36., -60., -36., 0.])) assert_array_equal(gy[gx.g.ic,:], np.array([0., 36., 60., 36., 12., -12., -36., -60., -36., 0.]))
def preevolve(self): """ preevolve is called before we being the timestepping loop. For the incompressible solver, this does an initial projection on the velocity field and then goes through the full evolution to get the value of phi. The fluid state (u, v) is then reset to values before this evolve. """ self.in_preevolve = True myg = self.cc_data.grid u = self.cc_data.get_var("x-velocity") v = self.cc_data.get_var("y-velocity") self.cc_data.fill_BC("x-velocity") self.cc_data.fill_BC("y-velocity") # 1. do the initial projection. This makes sure that our original # velocity field satisties div U = 0 # next create the multigrid object. We want Neumann BCs on phi # at solid walls and periodic on phi for periodic BCs mg = MG.CellCenterMG2d(myg.nx, myg.ny, xl_BC_type="periodic", xr_BC_type="periodic", yl_BC_type="periodic", yr_BC_type="periodic", xmin=myg.xmin, xmax=myg.xmax, ymin=myg.ymin, ymax=myg.ymax, verbose=0) # first compute divU divU = mg.soln_grid.scratch_array() divU.v()[:, :] = \ 0.5*(u.ip(1) - u.ip(-1))/myg.dx + 0.5*(v.jp(1) - v.jp(-1))/myg.dy # solve L phi = DU # initialize our guess to the solution, set the RHS to divU and # solve mg.init_zeros() mg.init_RHS(divU) mg.solve(rtol=1.e-10) # store the solution in our self.cc_data object -- include a single # ghostcell phi = self.cc_data.get_var("phi") phi[:, :] = mg.get_solution(grid=myg) # compute the cell-centered gradient of phi and update the # velocities gradp_x, gradp_y = mg.get_solution_gradient(grid=myg) u[:, :] -= gradp_x v[:, :] -= gradp_y # fill the ghostcells self.cc_data.fill_BC("x-velocity") self.cc_data.fill_BC("y-velocity") # 2. now get an approximation to gradp at n-1/2 by going through the # evolution. # store the current solution -- we'll restore it in a bit orig_data = patch.cell_center_data_clone(self.cc_data) # get the timestep self.method_compute_timestep() # evolve self.evolve() # update gradp_x and gradp_y in our main data object new_gp_x = self.cc_data.get_var("gradp_x") new_gp_y = self.cc_data.get_var("gradp_y") orig_gp_x = orig_data.get_var("gradp_x") orig_gp_y = orig_data.get_var("gradp_y") orig_gp_x[:, :] = new_gp_x[:, :] orig_gp_y[:, :] = new_gp_y[:, :] self.cc_data = orig_data if self.verbose > 0: print("done with the pre-evolution") self.in_preevolve = False
def evolve(self): """ Evolve the incompressible equations through one timestep. """ u = self.cc_data.get_var("x-velocity") v = self.cc_data.get_var("y-velocity") gradp_x = self.cc_data.get_var("gradp_x") gradp_y = self.cc_data.get_var("gradp_y") phi = self.cc_data.get_var("phi") myg = self.cc_data.grid # --------------------------------------------------------------------- # create the limited slopes of u and v (in both directions) # --------------------------------------------------------------------- limiter = self.rp.get_param("incompressible.limiter") ldelta_ux = reconstruction.limit(u, myg, 1, limiter) ldelta_vx = reconstruction.limit(v, myg, 1, limiter) ldelta_uy = reconstruction.limit(u, myg, 2, limiter) ldelta_vy = reconstruction.limit(v, myg, 2, limiter) # --------------------------------------------------------------------- # get the advective velocities # --------------------------------------------------------------------- """ the advective velocities are the normal velocity through each cell interface, and are defined on the cell edges, in a MAC type staggered form n+1/2 v i,j+1/2 +------+------+ | | n+1/2 | | n+1/2 u + U + u i-1/2,j | i,j | i+1/2,j | | +------+------+ n+1/2 v i,j-1/2 """ # this returns u on x-interfaces and v on y-interfaces. These # constitute the MAC grid if self.verbose > 0: print(" making MAC velocities") _um, _vm = incomp_interface.mac_vels(myg.ng, myg.dx, myg.dy, self.dt, u, v, ldelta_ux, ldelta_vx, ldelta_uy, ldelta_vy, gradp_x, gradp_y) u_MAC = ai.ArrayIndexer(d=_um, grid=myg) v_MAC = ai.ArrayIndexer(d=_vm, grid=myg) # --------------------------------------------------------------------- # do a MAC projection ot make the advective velocities divergence # free # --------------------------------------------------------------------- # we will solve L phi = D U^MAC, where phi is cell centered, and # U^MAC is the MAC-type staggered grid of the advective # velocities. if self.verbose > 0: print(" MAC projection") # create the multigrid object mg = MG.CellCenterMG2d(myg.nx, myg.ny, xl_BC_type="periodic", xr_BC_type="periodic", yl_BC_type="periodic", yr_BC_type="periodic", xmin=myg.xmin, xmax=myg.xmax, ymin=myg.ymin, ymax=myg.ymax, verbose=0) # first compute divU divU = mg.soln_grid.scratch_array() # MAC velocities are edge-centered. divU is cell-centered. divU.v()[:, :] = \ (u_MAC.ip(1) - u_MAC.v())/myg.dx + (v_MAC.jp(1) - v_MAC.v())/myg.dy # solve the Poisson problem mg.init_zeros() mg.init_RHS(divU) mg.solve(rtol=1.e-12) # update the normal velocities with the pressure gradient -- these # constitute our advective velocities phi_MAC = self.cc_data.get_var("phi-MAC") solution = mg.get_solution() phi_MAC.v(buf=1)[:, :] = solution.v(buf=1) # we need the MAC velocities on all edges of the computational domain b = (0, 1, 0, 0) u_MAC.v( buf=b)[:, :] -= (phi_MAC.v(buf=b) - phi_MAC.ip(-1, buf=b)) / myg.dx b = (0, 0, 0, 1) v_MAC.v( buf=b)[:, :] -= (phi_MAC.v(buf=b) - phi_MAC.jp(-1, buf=b)) / myg.dy # --------------------------------------------------------------------- # recompute the interface states, using the advective velocity # from above # --------------------------------------------------------------------- if self.verbose > 0: print(" making u, v edge states") _ux, _vx, _uy, _vy = \ incomp_interface.states(myg.ng, myg.dx, myg.dy, self.dt, u, v, ldelta_ux, ldelta_vx, ldelta_uy, ldelta_vy, gradp_x, gradp_y, u_MAC, v_MAC) u_xint = ai.ArrayIndexer(d=_ux, grid=myg) v_xint = ai.ArrayIndexer(d=_vx, grid=myg) u_yint = ai.ArrayIndexer(d=_uy, grid=myg) v_yint = ai.ArrayIndexer(d=_vy, grid=myg) # --------------------------------------------------------------------- # update U to get the provisional velocity field # --------------------------------------------------------------------- if self.verbose > 0: print(" doing provisional update of u, v") # compute (U.grad)U # we want u_MAC U_x + v_MAC U_y advect_x = myg.scratch_array() advect_y = myg.scratch_array() # u u_x + v u_y advect_x.v()[:, :] = \ 0.5*(u_MAC.v() + u_MAC.ip(1))*(u_xint.ip(1) - u_xint.v())/myg.dx + \ 0.5*(v_MAC.v() + v_MAC.jp(1))*(u_yint.jp(1) - u_yint.v())/myg.dy # u v_x + v v_y advect_y.v()[:, :] = \ 0.5*(u_MAC.v() + u_MAC.ip(1))*(v_xint.ip(1) - v_xint.v())/myg.dx + \ 0.5*(v_MAC.v() + v_MAC.jp(1))*(v_yint.jp(1) - v_yint.v())/myg.dy proj_type = self.rp.get_param("incompressible.proj_type") if proj_type == 1: u[:, :] -= (self.dt * advect_x[:, :] + self.dt * gradp_x[:, :]) v[:, :] -= (self.dt * advect_y[:, :] + self.dt * gradp_y[:, :]) elif proj_type == 2: u[:, :] -= self.dt * advect_x[:, :] v[:, :] -= self.dt * advect_y[:, :] self.cc_data.fill_BC("x-velocity") self.cc_data.fill_BC("y-velocity") # --------------------------------------------------------------------- # project the final velocity # --------------------------------------------------------------------- # now we solve L phi = D (U* /dt) if self.verbose > 0: print(" final projection") # create the multigrid object mg = MG.CellCenterMG2d(myg.nx, myg.ny, xl_BC_type="periodic", xr_BC_type="periodic", yl_BC_type="periodic", yr_BC_type="periodic", xmin=myg.xmin, xmax=myg.xmax, ymin=myg.ymin, ymax=myg.ymax, verbose=0) # first compute divU # u/v are cell-centered, divU is cell-centered divU.v()[:, :] = \ 0.5*(u.ip(1) - u.ip(-1))/myg.dx + 0.5*(v.jp(1) - v.jp(-1))/myg.dy mg.init_RHS(divU / self.dt) # use the old phi as our initial guess phiGuess = mg.soln_grid.scratch_array() phiGuess.v(buf=1)[:, :] = phi.v(buf=1) mg.init_solution(phiGuess) # solve mg.solve(rtol=1.e-12) # store the solution phi[:, :] = mg.get_solution(grid=myg) # compute the cell-centered gradient of p and update the velocities # this differs depending on what we projected. gradphi_x, gradphi_y = mg.get_solution_gradient(grid=myg) # u = u - grad_x phi dt u[:, :] -= self.dt * gradphi_x v[:, :] -= self.dt * gradphi_y # store gradp for the next step if proj_type == 1: gradp_x[:, :] += gradphi_x[:, :] gradp_y[:, :] += gradphi_y[:, :] elif proj_type == 2: gradp_x[:, :] = gradphi_x[:, :] gradp_y[:, :] = gradphi_y[:, :] self.cc_data.fill_BC("x-velocity") self.cc_data.fill_BC("y-velocity") if self.particles is not None: self.particles.update_particles(self.dt) # increment the time if not self.in_preevolve: self.cc_data.t += self.dt self.n += 1
def test_poisson_dirichlet(N, store_bench=False, comp_bench=False, make_plot=False, verbose=1): # test the multigrid solver nx = N ny = nx # create the multigrid object a = MG.CellCenterMG2d(nx, ny, xl_BC_type="dirichlet", yl_BC_type="dirichlet", xr_BC_type="dirichlet", yr_BC_type="dirichlet", verbose=verbose) # initialize the solution to 0 a.init_zeros() # initialize the RHS using the function f rhs = f(a.x2d, a.y2d) a.init_RHS(rhs) # solve to a relative tolerance of 1.e-11 a.solve(rtol=1.e-11) # alternately, we can just use smoothing by uncommenting the following #a.smooth(a.nlevels-1,50000) # get the solution v = a.get_solution() # compute the error from the analytic solution b = true(a.x2d, a.y2d) e = v - b print(" L2 error from true solution = %g\n rel. err from previous cycle = %g\n num. cycles = %d" % \ (e.norm(), a.relative_error, a.num_cycles)) # plot it if make_plot: plt.figure(num=1, figsize=(5.0, 5.0), dpi=100, facecolor='w') plt.imshow(np.transpose(v[a.ilo:a.ihi + 1, a.jlo:a.jhi + 1]), interpolation="nearest", origin="lower", extent=[a.xmin, a.xmax, a.ymin, a.ymax]) plt.xlabel("x") plt.ylabel("y") plt.savefig("mg_test.png") # store the output for later comparison bench = "mg_poisson_dirichlet" bench_dir = os.environ["PYRO_HOME"] + "/multigrid/tests/" my_data = a.get_solution_object() if store_bench: my_data.write("{}/{}".format(bench_dir, bench)) # do we do a comparison? if comp_bench: compare_file = "{}/{}".format(bench_dir, bench) msg.warning("comparing to: %s " % (compare_file)) bench_grid, bench_data = patch.read(compare_file) result = compare.compare(my_data.grid, my_data, bench_grid, bench_data) if result == 0: msg.success("results match benchmark\n") else: msg.warning("ERROR: " + compare.errors[result] + "\n") return result return None
def doit(nx, ny): nproj = 2 # create a mesh containing the x- and y-velocities, and periodic boundary # conditions myg = patch.Grid2d(nx, ny, ng=1) bc = bnd.BC(xlb="periodic", xrb="periodic", ylb="periodic", yrb="periodic") U = patch.CellCenterData2d(myg) U.register_var('u-old', bc) U.register_var('v-old', bc) U.register_var('u+gphi', bc) U.register_var('v+gphi', bc) U.register_var('u', bc) U.register_var('v', bc) U.register_var('divU', bc) U.register_var('phi-old', bc) U.register_var('phi', bc) U.register_var('dphi', bc) U.register_var('gradphi_x-old', bc) U.register_var('gradphi_y-old', bc) U.register_var('gradphi_x', bc) U.register_var('gradphi_y', bc) U.create() # initialize a divergence free velocity field, # u = -sin^2(pi x) sin(2 pi y), v = sin^2(pi y) sin(2 pi x) u = U.get_var('u') v = U.get_var('v') u[:, :] = -(np.sin(np.pi * myg.x2d)**2) * np.sin(2.0 * np.pi * myg.y2d) v[:, :] = (np.sin(np.pi * myg.y2d)**2) * np.sin(2.0 * np.pi * myg.x2d) # store the original, divergence free velocity field for comparison later uold = U.get_var('u-old') vold = U.get_var('v-old') uold[:, :] = u.copy() vold[:, :] = v.copy() # the projection routine should decompose U into a divergence free # part, U_d, plus the gradient of a scalar. Add on the gradient # of a scalar that satisfies gradphi.n = 0. After the projection, # we should recover the divergence free field above. Take phi to # be a gaussian, exp(-((x-x0)^2 + (y-y0)^2)/R) R = 0.1 x0 = 0.5 y0 = 0.5 phi = U.get_var('phi-old') gradphi_x = U.get_var('gradphi_x-old') gradphi_y = U.get_var('gradphi_y-old') phi[:, :] = np.exp(-((myg.x2d - x0)**2 + (myg.y2d - y0)**2) / R**2) gradphi_x[:, :] = phi * (-2.0 * (myg.x2d - x0) / R**2) gradphi_y[:, :] = phi * (-2.0 * (myg.y2d - y0) / R**2) u += gradphi_x v += gradphi_y u_plus_gradphi = U.get_var('u+gphi') v_plus_gradphi = U.get_var('v+gphi') u_plus_gradphi[:, :] = u[:, :] v_plus_gradphi[:, :] = v[:, :] # use the mesh class to enforce the periodic BCs on the velocity field U.fill_BC_all() # now compute the cell-centered, centered-difference divergence divU = U.get_var('divU') divU[myg.ilo:myg.ihi+1, myg.jlo:myg.jhi+1] = \ 0.5*(u[myg.ilo+1:myg.ihi+2, myg.jlo:myg.jhi+1] - u[myg.ilo-1:myg.ihi, myg.jlo:myg.jhi+1])/myg.dx + \ 0.5*(v[myg.ilo:myg.ihi+1, myg.jlo+1:myg.jhi+2] - v[myg.ilo:myg.ihi+1, myg.jlo-1:myg.jhi])/myg.dy # create the multigrid object with Neumann BCs a = MG.CellCenterMG2d(nx, ny, xl_BC_type="periodic", xr_BC_type="periodic", yl_BC_type="periodic", yr_BC_type="periodic", verbose=1) #-------------------------------------------------------------------------- # projections #-------------------------------------------------------------------------- for iproj in range(nproj): a.init_zeros() a.init_RHS(divU) a.solve(rtol=1.e-12) phi = U.get_var('phi') solution = a.get_solution() phi[myg.ilo-1:myg.ihi+2, myg.jlo-1:myg.jhi+2] = \ solution[a.ilo-1:a.ihi+2, a.jlo-1:a.jhi+2] dphi = U.get_var('dphi') dphi[:, :] = phi - U.get_var('phi-old') # compute the gradient of phi using centered differences gradphi_x = U.get_var('gradphi_x') gradphi_y = U.get_var('gradphi_y') gradphi_x[myg.ilo:myg.ihi+1, myg.jlo:myg.jhi+1] = \ 0.5*(phi[myg.ilo+1:myg.ihi+2, myg.jlo:myg.jhi+1] - phi[myg.ilo-1:myg.ihi, myg.jlo:myg.jhi+1])/myg.dx gradphi_y[myg.ilo:myg.ihi+1, myg.jlo:myg.jhi+1] = \ 0.5*(phi[myg.ilo:myg.ihi+1, myg.jlo+1:myg.jhi+2] - phi[myg.ilo:myg.ihi+1, myg.jlo-1:myg.jhi])/myg.dy # update the velocity field u -= gradphi_x v -= gradphi_y U.fill_BC_all() # recompute the divergence diagnostic divU[myg.ilo:myg.ihi+1, myg.jlo:myg.jhi+1] = \ 0.5*(u[myg.ilo+1:myg.ihi+2, myg.jlo:myg.jhi+1] - u[myg.ilo-1:myg.ihi, myg.jlo:myg.jhi+1])/myg.dx + \ 0.5*(v[myg.ilo:myg.ihi+1, myg.jlo+1:myg.jhi+2] - v[myg.ilo:myg.ihi+1, myg.jlo-1:myg.jhi])/myg.dy U.write("proj-periodic.after" + ("%d" % iproj))
def evolve(self, dt): """ Evolve the incompressible equations through one timestep. """ u = self.cc_data.get_var("x-velocity") v = self.cc_data.get_var("y-velocity") gradp_x = self.cc_data.get_var("gradp_x") gradp_y = self.cc_data.get_var("gradp_y") phi = self.cc_data.get_var("phi") myg = self.cc_data.grid dtdx = dt / myg.dx dtdy = dt / myg.dy #--------------------------------------------------------------------- # create the limited slopes of u and v (in both directions) #--------------------------------------------------------------------- limiter = self.rp.get_param("incompressible.limiter") if (limiter == 0): limitFunc = reconstruction_f.nolimit elif (limiter == 1): limitFunc = reconstruction_f.limit2 else: limitFunc = reconstruction_f.limit4 ldelta_ux = limitFunc(1, u, myg.qx, myg.qy, myg.ng) ldelta_vx = limitFunc(1, v, myg.qx, myg.qy, myg.ng) ldelta_uy = limitFunc(2, u, myg.qx, myg.qy, myg.ng) ldelta_vy = limitFunc(2, v, myg.qx, myg.qy, myg.ng) #--------------------------------------------------------------------- # get the advective velocities #--------------------------------------------------------------------- """ the advective velocities are the normal velocity through each cell interface, and are defined on the cell edges, in a MAC type staggered form n+1/2 v i,j+1/2 +------+------+ | | n+1/2 | | n+1/2 u + U + u i-1/2,j | i,j | i+1/2,j | | +------+------+ n+1/2 v i,j-1/2 """ # this returns u on x-interfaces and v on y-interfaces. These # constitute the MAC grid print(" making MAC velocities") u_MAC, v_MAC = incomp_interface_f.mac_vels(myg.qx, myg.qy, myg.ng, myg.dx, myg.dy, dt, u, v, ldelta_ux, ldelta_vx, ldelta_uy, ldelta_vy, gradp_x, gradp_y) #--------------------------------------------------------------------- # do a MAC projection ot make the advective velocities divergence # free #--------------------------------------------------------------------- # we will solve L phi = D U^MAC, where phi is cell centered, and # U^MAC is the MAC-type staggered grid of the advective # velocities. print(" MAC projection") # create the multigrid object mg = MG.CellCenterMG2d(myg.nx, myg.ny, xl_BC_type=self.cc_data.BCs["phi-MAC"].xlb, xr_BC_type=self.cc_data.BCs["phi-MAC"].xrb, yl_BC_type=self.cc_data.BCs["phi-MAC"].ylb, yr_BC_type=self.cc_data.BCs["phi-MAC"].yrb, xmin=myg.xmin, xmax=myg.xmax, ymin=myg.ymin, ymax=myg.ymax, verbose=0) # first compute divU divU = mg.soln_grid.scratch_array() # MAC velocities are edge-centered. divU is cell-centered. divU[mg.ilo:mg.ihi+1,mg.jlo:mg.jhi+1] = \ (u_MAC[myg.ilo+1:myg.ihi+2,myg.jlo:myg.jhi+1] - u_MAC[myg.ilo :myg.ihi+1,myg.jlo:myg.jhi+1])/myg.dx + \ (v_MAC[myg.ilo:myg.ihi+1,myg.jlo+1:myg.jhi+2] - v_MAC[myg.ilo:myg.ihi+1,myg.jlo :myg.jhi+1])/myg.dy # solve the Poisson problem mg.init_zeros() mg.init_RHS(divU) mg.solve(rtol=1.e-12) # update the normal velocities with the pressure gradient -- these # constitute our advective velocities phi_MAC = self.cc_data.get_var("phi-MAC") solution = mg.get_solution() phi_MAC[myg.ilo-1:myg.ihi+2,myg.jlo-1:myg.jhi+2] = \ solution[mg.ilo-1:mg.ihi+2,mg.jlo-1:mg.jhi+2] # we need the MAC velocities on all edges of the computational domain u_MAC[myg.ilo:myg.ihi+2,myg.jlo:myg.jhi+1] -= \ (phi_MAC[myg.ilo :myg.ihi+2,myg.jlo:myg.jhi+1] - phi_MAC[myg.ilo-1:myg.ihi+1,myg.jlo:myg.jhi+1])/myg.dx v_MAC[myg.ilo:myg.ihi+1,myg.jlo:myg.jhi+2] -= \ (phi_MAC[myg.ilo:myg.ihi+1,myg.jlo :myg.jhi+2] - phi_MAC[myg.ilo:myg.ihi+1,myg.jlo-1:myg.jhi+1])/myg.dy #--------------------------------------------------------------------- # recompute the interface states, using the advective velocity # from above #--------------------------------------------------------------------- print(" making u, v edge states") u_xint, v_xint, u_yint, v_yint = \ incomp_interface_f.states(myg.qx, myg.qy, myg.ng, myg.dx, myg.dy, dt, u, v, ldelta_ux, ldelta_vx, ldelta_uy, ldelta_vy, gradp_x, gradp_y, u_MAC, v_MAC) #--------------------------------------------------------------------- # update U to get the provisional velocity field #--------------------------------------------------------------------- print(" doing provisional update of u, v") # compute (U.grad)U # we want u_MAC U_x + v_MAC U_y advect_x = myg.scratch_array() advect_y = myg.scratch_array() advect_x[myg.ilo:myg.ihi+1,myg.jlo:myg.jhi+1] = \ 0.5*(u_MAC[myg.ilo :myg.ihi+1,myg.jlo:myg.jhi+1] + u_MAC[myg.ilo+1:myg.ihi+2,myg.jlo:myg.jhi+1]) * \ (u_xint[myg.ilo+1:myg.ihi+2,myg.jlo:myg.jhi+1] - u_xint[myg.ilo :myg.ihi+1,myg.jlo:myg.jhi+1])/myg.dx + \ 0.5*(v_MAC[myg.ilo:myg.ihi+1,myg.jlo :myg.jhi+1] + v_MAC[myg.ilo:myg.ihi+1,myg.jlo+1:myg.jhi+2]) * \ (u_yint[myg.ilo:myg.ihi+1,myg.jlo+1:myg.jhi+2] - u_yint[myg.ilo:myg.ihi+1,myg.jlo :myg.jhi+1])/myg.dy advect_y[myg.ilo:myg.ihi+1,myg.jlo:myg.jhi+1] = \ 0.5*(u_MAC[myg.ilo :myg.ihi+1,myg.jlo:myg.jhi+1] + u_MAC[myg.ilo+1:myg.ihi+2,myg.jlo:myg.jhi+1]) * \ (v_xint[myg.ilo+1:myg.ihi+2,myg.jlo:myg.jhi+1] - v_xint[myg.ilo :myg.ihi+1,myg.jlo:myg.jhi+1])/myg.dx + \ 0.5*(v_MAC[myg.ilo:myg.ihi+1,myg.jlo :myg.jhi+1] + v_MAC[myg.ilo:myg.ihi+1,myg.jlo+1:myg.jhi+2]) * \ (v_yint[myg.ilo:myg.ihi+1,myg.jlo+1:myg.jhi+2] - v_yint[myg.ilo:myg.ihi+1,myg.jlo :myg.jhi+1])/myg.dy proj_type = self.rp.get_param("incompressible.proj_type") if (proj_type == 1): u[:, :] -= (dt * advect_x[:, :] + dt * gradp_x[:, :]) v[:, :] -= (dt * advect_y[:, :] + dt * gradp_y[:, :]) elif (proj_type == 2): u[:, :] -= dt * advect_x[:, :] v[:, :] -= dt * advect_y[:, :] self.cc_data.fill_BC("x-velocity") self.cc_data.fill_BC("y-velocity") #--------------------------------------------------------------------- # project the final velocity #--------------------------------------------------------------------- # now we solve L phi = D (U* /dt) print(" final projection") # create the multigrid object mg = MG.CellCenterMG2d(myg.nx, myg.ny, xl_BC_type=self.cc_data.BCs["phi"].xlb, xr_BC_type=self.cc_data.BCs["phi"].xrb, yl_BC_type=self.cc_data.BCs["phi"].ylb, yr_BC_type=self.cc_data.BCs["phi"].yrb, xmin=myg.xmin, xmax=myg.xmax, ymin=myg.ymin, ymax=myg.ymax, verbose=0) # first compute divU # u/v are cell-centered, divU is cell-centered divU[mg.ilo:mg.ihi+1,mg.jlo:mg.jhi+1] = \ 0.5*(u[myg.ilo+1:myg.ihi+2,myg.jlo:myg.jhi+1] - u[myg.ilo-1:myg.ihi ,myg.jlo:myg.jhi+1])/myg.dx + \ 0.5*(v[myg.ilo:myg.ihi+1,myg.jlo+1:myg.jhi+2] - v[myg.ilo:myg.ihi+1,myg.jlo-1:myg.jhi ])/myg.dy mg.init_RHS(divU / dt) # use the old phi as our initial guess phiGuess = mg.soln_grid.scratch_array() phiGuess[mg.ilo-1:mg.ihi+2,mg.jlo-1:mg.jhi+2] = \ phi[myg.ilo-1:myg.ihi+2,myg.jlo-1:myg.jhi+2] mg.init_solution(phiGuess) # solve mg.solve(rtol=1.e-12) # store the solution solution = mg.get_solution() phi[myg.ilo-1:myg.ihi+2,myg.jlo-1:myg.jhi+2] = \ solution[mg.ilo-1:mg.ihi+2,mg.jlo-1:mg.jhi+2] # compute the cell-centered gradient of p and update the velocities # this differs depending on what we projected. gradphi_x = myg.scratch_array() gradphi_y = myg.scratch_array() gradphi_x[myg.ilo:myg.ihi+1,myg.jlo:myg.jhi+1] = \ 0.5*(phi[myg.ilo+1:myg.ihi+2,myg.jlo:myg.jhi+1] - phi[myg.ilo-1:myg.ihi ,myg.jlo:myg.jhi+1])/myg.dx gradphi_y[myg.ilo:myg.ihi+1,myg.jlo:myg.jhi+1] = \ 0.5*(phi[myg.ilo:myg.ihi+1,myg.jlo+1:myg.jhi+2] - phi[myg.ilo:myg.ihi+1,myg.jlo-1:myg.jhi ])/myg.dy # u = u - grad_x phi dt u[:, :] -= dt * gradphi_x v[:, :] -= dt * gradphi_y # store gradp for the next step if (proj_type == 1): gradp_x[:, :] += gradphi_x[:, :] gradp_y[:, :] += gradphi_y[:, :] elif (proj_type == 2): gradp_x[:, :] = gradphi_x[:, :] gradp_y[:, :] = gradphi_y[:, :] self.cc_data.fill_BC("x-velocity") self.cc_data.fill_BC("y-velocity")
def preevolve(self): """ preevolve is called before we being the timestepping loop. For the incompressible solver, this does an initial projection on the velocity field and then goes through the full evolution to get the value of phi. The fluid state (u, v) is then reset to values before this evolve. """ myg = self.cc_data.grid u = self.cc_data.get_var("x-velocity") v = self.cc_data.get_var("y-velocity") self.cc_data.fill_BC("x-velocity") self.cc_data.fill_BC("y-velocity") # 1. do the initial projection. This makes sure that our original # velocity field satisties div U = 0 # next create the multigrid object. We defined phi with # the right BCs previously mg = MG.CellCenterMG2d(myg.nx, myg.ny, xl_BC_type=self.cc_data.BCs["phi"].xlb, xr_BC_type=self.cc_data.BCs["phi"].xrb, yl_BC_type=self.cc_data.BCs["phi"].ylb, yr_BC_type=self.cc_data.BCs["phi"].yrb, xmin=myg.xmin, xmax=myg.xmax, ymin=myg.ymin, ymax=myg.ymax, verbose=0) # first compute divU divU = mg.soln_grid.scratch_array() divU[mg.ilo:mg.ihi+1,mg.jlo:mg.jhi+1] = \ 0.5*(u[myg.ilo+1:myg.ihi+2,myg.jlo:myg.jhi+1] - u[myg.ilo-1:myg.ihi ,myg.jlo:myg.jhi+1])/myg.dx + \ 0.5*(v[myg.ilo:myg.ihi+1,myg.jlo+1:myg.jhi+2] - v[myg.ilo:myg.ihi+1,myg.jlo-1:myg.jhi ])/myg.dy # solve L phi = DU # initialize our guess to the solution, set the RHS to divU and # solve mg.init_zeros() mg.init_RHS(divU) mg.solve(rtol=1.e-10) # store the solution in our self.cc_data object -- include a single # ghostcell phi = self.cc_data.get_var("phi") solution = mg.get_solution() phi[myg.ilo-1:myg.ihi+2,myg.jlo-1:myg.jhi+2] = \ solution[mg.ilo-1:mg.ihi+2,mg.jlo-1:mg.jhi+2] # compute the cell-centered gradient of phi and update the # velocities gradp_x = myg.scratch_array() gradp_y = myg.scratch_array() gradp_x[myg.ilo:myg.ihi+1,myg.jlo:myg.jhi+1] = \ 0.5*(phi[myg.ilo+1:myg.ihi+2,myg.jlo:myg.jhi+1] - phi[myg.ilo-1:myg.ihi ,myg.jlo:myg.jhi+1])/myg.dx gradp_y[myg.ilo:myg.ihi+1,myg.jlo:myg.jhi+1] = \ 0.5*(phi[myg.ilo:myg.ihi+1,myg.jlo+1:myg.jhi+2] - phi[myg.ilo:myg.ihi+1,myg.jlo-1:myg.jhi ])/myg.dy u[:, :] -= gradp_x v[:, :] -= gradp_y # fill the ghostcells self.cc_data.fill_BC("x-velocity") self.cc_data.fill_BC("y-velocity") # 2. now get an approximation to gradp at n-1/2 by going through the # evolution. # store the current solution -- we'll restore it in a bit orig_data = patch.cell_center_data_clone(self.cc_data) # get the timestep dt = self.timestep() # evolve self.evolve(dt) # update gradp_x and gradp_y in our main data object new_gp_x = self.cc_data.get_var("gradp_x") new_gp_y = self.cc_data.get_var("gradp_y") orig_gp_x = orig_data.get_var("gradp_x") orig_gp_y = orig_data.get_var("gradp_y") orig_gp_x[:, :] = new_gp_x[:, :] orig_gp_y[:, :] = new_gp_y[:, :] self.cc_data = orig_data print("done with the pre-evolution")