def crossflow_plume(fig): """ Define, run, and plot the simulations for a pure bubble plume in crossflow for validation to data in Socolofsky and Adams (2002). """ # Jet initial conditions z0 = 0.64 U0 = 0. phi_0 = -np.pi / 2. theta_0 = 0. D = 0.01 Tj = 21. + 273.15 Sj = 0. cj = 1. chem_name = 'tracer' # Ambient conditions ua = 0.15 T = 0. F = 0. H = 1.0 # Create the correct ambient profile data uj = U0 * np.cos(phi_0) * np.cos(theta_0) vj = U0 * np.cos(phi_0) * np.sin(theta_0) wj = U0 * np.sin(phi_0) profile_fname = './crossflow_plume.nc' profile = get_profile(profile_fname, z0, D, uj, vj, wj, Tj, Sj, ua, T, F, 1., H) # Create a bent plume model simulation object jlm = bpm.Model(profile) # Define the dispersed phase input to the model composition = ['nitrogen', 'oxygen', 'argon', 'carbon_dioxide'] mol_frac = np.array([0.78084, 0.20946, 0.009340, 0.00036]) air = dbm.FluidParticle(composition) particles = [] # Large bubbles Q_N = 0.5 / 60. / 1000. de0 = 0.008 T0 = Tj lambda_1 = 1. (m0, T0, nb0, P, Sa, Ta) = dispersed_phases.initial_conditions(profile, z0, air, mol_frac, Q_N, 1, de0, T0) particles.append( bpm.Particle(0., 0., z0, air, m0, T0, nb0, lambda_1, P, Sa, Ta, K=1., K_T=1., fdis=1.e-6)) # Small bubbles Q_N = 0.5 / 60. / 1000. de0 = 0.0003 T0 = Tj lambda_1 = 1. (m0, T0, nb0, P, Sa, Ta) = dispersed_phases.initial_conditions(profile, z0, air, mol_frac, Q_N, 1, de0, T0) particles.append( bpm.Particle(0., 0., z0, air, m0, T0, nb0, lambda_1, P, Sa, Ta, K=1., K_T=1., fdis=1.e-6)) # Run the simulation jlm.simulate(np.array([0., 0., z0]), D, U0, phi_0, theta_0, Sj, Tj, cj, chem_name, particles, track=True, dt_max=60., sd_max=100.) # Perpare variables for plotting xp = jlm.q[:, 7] / jlm.D yp = jlm.q[:, 9] / jlm.D plt.figure(fig) plt.clf() plt.show() ax1 = plt.subplot(111) ax1.plot(xp, yp, 'b-') ax1.set_xlabel('x / D') ax1.set_ylabel('z / D') ax1.invert_yaxis() ax1.grid(b=True, which='major', color='0.65', linestyle='-') plt.draw() return jlm
def get_particles(self, composition, data, md_gas0, md_oil0, profile, d50_gas, d50_oil, nbins, T0, z0, dispersant, sigma_fac, oil, mass_frac, hydrate, inert_drop): """ docstring for get_particles """ # Reduce surface tension if dispersant is applied if dispersant is True: sigma = np.array([[1.], [1.]]) * sigma_fac else: sigma = np.array([[1.], [1.]]) # Create DBM objects for the bubbles and droplets bubl = dbm.FluidParticle(composition, fp_type=0, sigma_correction=sigma[0], user_data=data) drop = dbm.FluidParticle(composition, fp_type=1, sigma_correction=sigma[1], user_data=data) # Get the local ocean conditions T, S, P = profile.get_values(z0, ['temperature', 'salinity', 'pressure']) rho = seawater.density(T, S, P) # Get the mole fractions of the released fluids molf_gas = bubl.mol_frac(md_gas0) molf_oil = drop.mol_frac(md_oil0) print molf_gas print molf_oil # Use the Rosin-Rammler distribution to get the mass flux in each # size class # de_gas, md_gas = sintef.rosin_rammler(nbins, d50_gas, np.sum(md_gas0), # bubl.interface_tension(md_gas0, T0, S, P), # bubl.density(md_gas0, T0, P), rho) # de_oil, md_oil = sintef.rosin_rammler(nbins, d50_oil, np.sum(md_oil0), # drop.interface_tension(md_oil0, T0, S, P), # drop.density(md_oil0, T0, P), rho) # Get the user defined particle size distibution de_oil, vf_oil, de_gas, vf_gas = self.userdefined_de() md_gas = np.sum(md_gas0) * vf_gas md_oil = np.sum(md_oil0) * vf_oil # Define a inert particle to be used if inert liquid particles are use # in the simulations molf_inert = 1. isfluid = True iscompressible = True rho_o = drop.density(md_oil0, T0, P) inert = dbm.InsolubleParticle(isfluid, iscompressible, rho_p=rho_o, gamma=40., beta=0.0007, co=2.90075e-9) # Create the particle objects particles = [] t_hyd = 0. # Bubbles for i in range(nbins): if md_gas[i] > 0.: (m0, T0, nb0, P, Sa, Ta) = dispersed_phases.initial_conditions( profile, z0, bubl, molf_gas, md_gas[i], 2, de_gas[i], T0) # Get the hydrate formation time for bubbles if hydrate is True and dispersant is False: t_hyd = dispersed_phases.hydrate_formation_time(bubl, z0, m0, T0, profile) if np.isinf(t_hyd): t_hyd = 0. else: t_hyd = 0. particles.append(bpm.Particle(0., 0., z0, bubl, m0, T0, nb0, 1.0, P, Sa, Ta, K=1., K_T=1., fdis=1.e-6, t_hyd=t_hyd)) # Droplets for i in range(len(de_oil)): # Add the live droplets to the particle list if md_oil[i] > 0. and not inert_drop: (m0, T0, nb0, P, Sa, Ta) = dispersed_phases.initial_conditions( profile, z0, drop, molf_oil, md_oil[i], 2, de_oil[i], T0) # Get the hydrate formation time for bubbles if hydrate is True and dispersant is False: t_hyd = dispersed_phases.hydrate_formation_time(drop, z0, m0, T0, profile) if np.isinf(t_hyd): t_hyd = 0. else: t_hyd = 0. particles.append(bpm.Particle(0., 0., z0, drop, m0, T0, nb0, 1.0, P, Sa, Ta, K=1., K_T=1., fdis=1.e-6, t_hyd=t_hyd)) # Add the inert droplets to the particle list if md_oil[i] > 0. and inert_drop is True: (m0, T0, nb0, P, Sa, Ta) = dispersed_phases.initial_conditions( profile, z0, inert, molf_oil, md_oil[i], 2, de_oil[i], T0) particles.append(bpm.Particle(0., 0., z0, inert, m0, T0, nb0, 1.0, P, Sa, Ta, K=1., K_T=1., fdis=1.e-6, t_hyd=0.)) # Define the lambda for particles model = params.Scales(profile, particles) for j in range(len(particles)): particles[j].lambda_1 = model.lambda_1(z0, j) # Return the particle list return particles
# Larger free gas bubbles mb0 = 5. # total mass flux in kg/s de = 0.005 # bubble diameter in m lambda_1 = 0.85 (m0, T0, nb0, P, Sa, Ta) = dispersed_phases.initial_conditions(ctd, z0, gas, yk, mb0, 2, de, Tj) disp_phases.append( bent_plume_model.Particle(0., 0., z0, gas, m0, T0, nb0, lambda_1, P, Sa, Ta, K=1., K_T=1., fdis=1.e-6, t_hyd=0., lag_time=False)) # Smaller free gas bubbles mb0 = 5. # total mass flux in kg/s de = 0.0005 # bubble diameter in m lambda_1 = 0.95 (m0, T0, nb0, P, Sa, Ta) = dispersed_phases.initial_conditions(ctd, z0, gas, yk, mb0, 2, de, Tj)
def particles(m_tot, d, vf, profile, oil, yk, x0, y0, z0, Tj, lambda_1, lag_time): """ Create particles to add to a bent plume model simulation Creates bent_plume_model.Particle objects for the given particle properties so that they can be added to the total list of particles in the simulation. Parameters ---------- m_tot : float Total mass flux of this fluid phase in the simulation (kg/s) d : np.array Array of particle sizes for this fluid phase (m) vf : np.array Array of volume fractions for each particle size for this fluid phase (--). This array should sum to 1.0. profile : ambient.Profile An ambient.Profile object with the ambient ocean water column data oil : dbm.FluidParticle A dbm.FluidParticle object that contains the desired oil database composition yk : np.array Mole fractions of each compound in the chemical database of the oil dbm.FluidParticle object (--). x0, y0, z0 : floats Initial position of the particles in the simulation domain (m). Note that x0 and y0 should be zero for particles starting on the plume centerline. Tj : float Initial temperature of the particles in the jet (K) lambda_1 : float Value of the dispersed phase spreading parameter of the jet integral model (--). lag_time : bool Flag that indicates whether (True) or not (False) to use the biodegradation lag times data. Returns ------- disp_phases : list of bent_plume_model.Particle objects List of `bent_plume_model.Particle` objects to be added to the present bent plume model simulation based on the given input data. Notes ----- See the documentation for the `bent_plume_model` for more information on the `Particle` object. """ # Create an empty list of particles disp_phases = [] # Add each particle in the distribution separately for i in range(len(d)): # Get the total mass flux of this fluid phase for the present # particle size mb0 = vf[i] * m_tot # Get the properties of these particles at the source (m0, T0, nb0, P, Sa, Ta) = dispersed_phases.initial_conditions(profile, z0, oil, yk, mb0, 2, d[i], Tj) # Append these particles to the list of particles in the simulation disp_phases.append( bent_plume_model.Particle(x0, y0, z0, oil, m0, T0, nb0, lambda_1, P, Sa, Ta, K=1., K_T=1., fdis=1.e-6, t_hyd=0., lag_time=lag_time)) # Return the list of particles return disp_phases
def get_sim_data(): """ Create the data needed to initialize a simulation Performs the steps necessary to set up a bent plume model simulation and passes the input variables to the `Model` object and `Model.simulate()` method. Returns ------- profile : `ambient.Profile` object Return a profile object from the BM54 CTD data z0 : float Depth of the release port (m) D : float Diameter of the release port (m) Vj : float Initial velocity of the jet (m/s) phi_0 : float Vertical angle from the horizontal for the discharge orientation (rad in range +/- pi/2) theta_0 : float Horizontal angle from the x-axis for the discharge orientation. The x-axis is taken in the direction of the ambient current. (rad in range 0 to 2 pi) Sj : float Salinity of the continuous phase fluid in the discharge (psu) Tj : float Temperature of the continuous phase fluid in the discharge (T) cj : ndarray Concentration of passive tracers in the discharge (user-defined) tracers : string list List of passive tracers in the discharge. These can be chemicals present in the ambient `profile` data, and if so, entrainment of these chemicals will change the concentrations computed for these tracers. However, none of these concentrations are used in the dissolution of the dispersed phase. Hence, `tracers` should not contain any chemicals present in the dispersed phase particles. particles : list of `Particle` objects List of `Particle` objects describing each dispersed phase in the simulation dt_max : float Maximum step size to take in the storage of the simulation solution (s) sd_max : float Maximum number of orifice diameters to compute the solution along the plume centerline (m/m) """ # Get the ambient CTD data profile = get_profile() # Specify the release location and geometry and initialize a particle # list z0 = 300. D = 0.3 particles = [] # Add a dissolving particle to the list composition = ['oxygen', 'nitrogen', 'argon'] yk = np.array([1.0, 0., 0.]) o2 = dbm.FluidParticle(composition) Q_N = 1.5 / 60. / 60. de = 0.009 lambda_1 = 0.85 (m0, T0, nb0, P, Sa, Ta) = dispersed_phases.initial_conditions( profile, z0, o2, yk, Q_N, 1, de) particles.append(bent_plume_model.Particle(0., 0., z0, o2, m0, T0, nb0, lambda_1, P, Sa, Ta, K=1., K_T=1., fdis=1.e-6, t_hyd=0.)) # Add an insoluble particle to the list composition = ['inert'] yk = np.array([1.]) oil = dbm.InsolubleParticle(True, True) mb0 = 1. de = 0.01 lambda_1 = 0.8 (m0, T0, nb0, P, Sa, Ta) = dispersed_phases.initial_conditions( profile, z0, oil, yk, mb0, 1, de) particles.append(bent_plume_model.Particle(0., 0., z0, oil, m0, T0, nb0, lambda_1, P, Sa, Ta, K=1., K_T=1., fdis=1.e-6, t_hyd=0.)) # Set the other simulation parameters Vj = 0. phi_0 = -np.pi/2. theta_0 = 0. Sj = 0. Tj = Ta cj = np.array([1.]) tracers = ['tracer'] dt_max = 60. sd_max = 3000. # Return the results return (profile, np.array([0., 0., z0]), D, Vj, phi_0, theta_0, Sj, Tj, cj, tracers, particles, dt_max, sd_max)
def test_particle_obj(): """ Test the object behavior for the `Particle` object Test the instantiation and attribute data for the `Particle` object of the `bent_plume_model` module. """ # Set up the base parameters describing a particle object T = 273.15 + 15. P = 150e5 Sa = 35. Ta = 273.15 + 4. composition = ['methane', 'ethane', 'propane', 'oxygen'] yk = np.array([0.85, 0.07, 0.08, 0.0]) de = 0.005 lambda_1 = 0.85 K = 1. Kt = 1. fdis = 1e-6 x = 0. y = 0. z = 0. # Compute a few derived quantities bub = dbm.FluidParticle(composition) nb0 = 1.e5 m0 = bub.masses_by_diameter(de, T, P, yk) # Create a `PlumeParticle` object bub_obj = bent_plume_model.Particle(x, y, z, bub, m0, T, nb0, lambda_1, P, Sa, Ta, K, Kt, fdis) # Check if the initialized attributes are correct assert bub_obj.integrate == True assert bub_obj.sim_stored == False assert bub_obj.farfield == False assert bub_obj.t == 0. assert bub_obj.x == x assert bub_obj.y == y assert bub_obj.z == z for i in range(len(composition)): assert bub_obj.composition[i] == composition[i] assert_array_almost_equal(bub_obj.m0, m0, decimal=6) assert bub_obj.T0 == T assert_array_almost_equal(bub_obj.m, m0, decimal=6) assert bub_obj.T == T assert bub_obj.cp == seawater.cp() * 0.5 assert bub_obj.K == K assert bub_obj.K_T == Kt assert bub_obj.fdis == fdis for i in range(len(composition)-1): assert bub_obj.diss_indices[i] == True assert bub_obj.diss_indices[-1] == False assert bub_obj.nb0 == nb0 assert bub_obj.lambda_1 == lambda_1 # Including the values after the first call to the update method us_ans = bub.slip_velocity(m0, T, P, Sa, Ta) rho_p_ans = bub.density(m0, T, P) A_ans = bub.surface_area(m0, T, P, Sa, Ta) Cs_ans = bub.solubility(m0, T, P, Sa) beta_ans = bub.mass_transfer(m0, T, P, Sa, Ta) beta_T_ans = bub.heat_transfer(m0, T, P, Sa, Ta) assert bub_obj.us == us_ans assert bub_obj.rho_p == rho_p_ans assert bub_obj.A == A_ans assert_array_almost_equal(bub_obj.Cs, Cs_ans, decimal=6) assert_array_almost_equal(bub_obj.beta, beta_ans, decimal=6) assert bub_obj.beta_T == beta_T_ans # Test the bub_obj.outside() method bub_obj.outside(Ta, Sa, P) assert bub_obj.us == 0. assert bub_obj.rho_p == seawater.density(Ta, Sa, P) assert bub_obj.A == 0. assert_array_almost_equal(bub_obj.Cs, np.zeros(len(composition))) assert_array_almost_equal(bub_obj.beta, np.zeros(len(composition))) assert bub_obj.beta_T == 0. assert bub_obj.T == Ta # No need to test the properties or diameter objects since they are # inherited from the `single_bubble_model` and tested in `test_sbm`. # No need to test the bub_obj.track(), bub_obj.run_sbm() since they will # be tested below for the simulation cases. # Check functionality of insoluble particle drop = dbm.InsolubleParticle(isfluid=True, iscompressible=True) m0 = drop.mass_by_diameter(de, T, P, Sa, Ta) drop_obj = bent_plume_model.Particle(x, y, z, drop, m0, T, nb0, lambda_1, P, Sa, Ta, K, fdis=fdis, K_T=Kt) assert len(drop_obj.composition) == 1 assert drop_obj.composition[0] == 'inert' assert_array_almost_equal(drop_obj.m0, m0, decimal=6) assert drop_obj.T0 == T assert_array_almost_equal(drop_obj.m, m0, decimal=6) assert drop_obj.T == T assert drop_obj.cp == seawater.cp() * 0.5 assert drop_obj.K == K assert drop_obj.K_T == Kt assert drop_obj.fdis == fdis assert drop_obj.diss_indices[0] == True assert drop_obj.nb0 == nb0 assert drop_obj.lambda_1 == lambda_1 # Including the values after the first call to the update method us_ans = drop.slip_velocity(m0, T, P, Sa, Ta) rho_p_ans = drop.density(T, P, Sa, Ta) A_ans = drop.surface_area(m0, T, P, Sa, Ta) beta_T_ans = drop.heat_transfer(m0, T, P, Sa, Ta) assert drop_obj.us == us_ans assert drop_obj.rho_p == rho_p_ans assert drop_obj.A == A_ans assert drop_obj.beta_T == beta_T_ans