/* * Stretch and tilt vortons using velocity field * timeStep - amount of time by which to advance simulation * uFrame - frame counter * * see J. T. Beale, A convergent three-dimensional vortex method with * grid-free stretching, Math. Comp. 46 (1986), 401-24, April. * * This routine assumes CreateInfluenceTree has already executed. * */ void StretchAndTiltVortons(float timeStep) { if ((0.0f == mVelGrid.GetExtent().x) || (0.0f == mVelGrid.GetExtent().y) || (0.0f == mVelGrid.GetExtent().z)) { // Domain is 2D, so stretching & tilting does not occur. return; } // Compute all gradients of all components of velocity. UniformGrid <Matrix3x3> velocityJacobianGrid = new UniformGrid <Matrix3x3>(mVelGrid); velocityJacobianGrid.Init(); UniformGridMath.ComputeJacobian(ref velocityJacobianGrid, mVelGrid); int numVortons = mVortons.Count; for (int offset = 0; offset < numVortons; ++offset) { // For each vorton... Matrix3x3 velJac = (Matrix3x3)velocityJacobianGrid.Interpolate(mVortons[offset].position); Vector3 stretchTilt = mVortons[offset].vorticity * velJac; // Usual way to compute stretching & tilting mVortons[offset].vorticity += /* fudge factor for stability */ 0.5f * stretchTilt * timeStep; } }
/* * Create base layer of vorton influence tree. * * This is the leaf layer, where each grid cell corresponds(on average) to * a single vorton.Some cells might contain multiple vortons and some zero. * * Each cell effectively has a single "supervorton" which its parent layers * in the influence tree will in turn aggregate. * * This implementation of gridifying the base layer is NOT suitable * for Eulerian operations like approximating spatial derivatives * * of vorticity or solving a vector Poisson equation, because this * routine associates each vortex with a single corner point of the * grid cell that contains it. To create a grid for Eulerian calculations, * each vorton would contribute to all 8 corner points of the grid * cell that contains it. * * We could rewrite this to suit "Eulerian" operations, in which case * we would want to omit "size" and "position" since the grid would * * implicitly represent that information. That concern goes hand-in-hand * with the method used to compute velocity from vorticity. * * Ultimately we need to make sure theoretically conserved quantities behave as expected. * * This method assumes the influence tree skeleton has already been created, * and the leaf layer initialized to all "zeros", meaning it contains no vortons. */ public void MakeBaseVortonGrid() { int numVortons = mVortons.Count; UniformGrid <VortonClusterAux> ugAux = new UniformGrid <VortonClusterAux>(mInfluenceTree[0]); // Temporary auxilliary information used during aggregation. ugAux.Init(); // Compute preliminary vorticity grid. for (int uVorton = 0; uVorton < numVortons; ++uVorton) { // For each vorton in this simulation... Vector3 position = mVortons[uVorton].position; uint uOffset = mInfluenceTree[0].OffsetOfPosition(position); float vortMag = mVortons[uVorton].vorticity.magnitude; mInfluenceTree[0][uOffset].position += mVortons[uVorton].position * vortMag; // Compute weighted position -- to be normalized later. mInfluenceTree[0][uOffset].vorticity += mVortons[uVorton].vorticity; // Tally vorticity sum. mInfluenceTree[0][uOffset].radius = mVortons[uVorton].radius; // Assign volume element size. // OBSOLETE. See comments below: UpdateBoundingBox( rVortonAux.mMinCorner , rVortonAux.mMaxCorner , rVorton.mPosition ) ; ugAux[uOffset].mVortNormSum += vortMag; // Accumulate vorticity on the VortonClusterAux } // Post-process preliminary grid (VortonClusterAux); normalize center-of-vorticity and compute sizes, for each grid cell. uint[] num = { mInfluenceTree[0].GetNumPoints(0), mInfluenceTree[0].GetNumPoints(1), mInfluenceTree[0].GetNumPoints(2) }; uint numXY = num[0] * num[1]; uint[] idx = new uint[3]; for (idx[2] = 0; idx[2] < num[2]; ++idx[2]) { uint zShift = idx[2] * numXY; for (idx[1] = 0; idx[1] < num[1]; ++idx[1]) { uint yzShift = idx[1] * num[0] + zShift; for (idx[0] = 0; idx[0] < num[0]; ++idx[0]) { uint offset = idx[0] + yzShift; VortonClusterAux rVortonAux = ugAux[offset]; if (rVortonAux.mVortNormSum != float.Epsilon) { // This cell contains at least one vorton. // Normalize weighted position sum to obtain center-of-vorticity. mInfluenceTree[0][offset].position /= rVortonAux.mVortNormSum; } } } } }
/* * Compute velocity due to vortons, for every point in a uniform grid * This routine assumes CreateInfluenceTree has already executed. */ void ComputeVelocityGrid() { mVelGrid.Clear(); // Clear any stale velocity information mVelGrid.CopyShape(mInfluenceTree[0]); // Use same shape as base vorticity grid. (Note: could differ if you want.) mVelGrid.Init(); // Reserve memory for velocity grid. uint numZ = mVelGrid.GetNumPoints(2); if (mUseMultithreads) { // Estimate grain size based on size of problem and number of processors. int grainSize = (int)Mathf.Max(1, numZ / numberOfProcessors); List <ManualResetEvent> handles = new List <ManualResetEvent>(); for (var i = 0; i < numberOfProcessors; i++) { ManualResetEvent handle = new ManualResetEvent(false); handles.Add(handle); // Send the custom object to the threaded method. ThreadInfo threadInfo = new ThreadInfo(); threadInfo.begin = i * grainSize; threadInfo.end = (i + 1) * grainSize; threadInfo.timeStep = 0; threadInfo.vortonSim = this; threadInfo.handle = handle; WaitCallback callBack = new WaitCallback(ComputeVelocityGridSliceThreaded); Nyahoon.ThreadPool.QueueUserWorkItem(callBack, threadInfo); } WaitHandle.WaitAll(handles.ToArray()); } else { ComputeVelocityGridSlice(0, numZ); } }
/* * Diffuse vorticity using a particle strength exchange method. * * This routine partitions space into cells using the same grid * as the "base vorton" grid. Each vorton gets assigned to the * cell that contains it. Then, each vorton exchanges some * of its vorticity with its neighbors in adjacent cells. * * This routine makes some simplifying assumptions to speed execution: * * - Distance does not influence the amount of vorticity exchanged, * except in as much as only vortons within a certain region of * each other exchange vorticity. This amounts to saying our kernel, * eta, is a top-hat function. * * - Theoretically, if an adjacent cell contains no vortons * then this simulation should generate vorticity within * that cell, e.g. by creating a new vorton in the adjacent cell. * * - This simulation reduces the vorticity of each vorton, alleging * that this vorticity is dissipated analogously to how energy * dissipates at Kolmogorov microscales. This treatment is not * realistic but it retains qualitative characteristics that we * want, e.g. that the flow dissipates at a rate related to viscosity. * Dissipation in real flows is a more complicated phenomenon. * * see Degond & Mas-Gallic (1989): The weighted particle method for * convection-diffusion equations, part 1: the case of an isotropic viscosity. * Math. Comput., v. 53, n. 188, pp. 485-507, October. * * timeStep - amount of time by which to advance simulation * * This routine assumes CreateInfluenceTree has already executed. * */ void DiffuseVorticityPSE(float timeStep) { // Phase 1: Partition vortons // Create a spatial partition for the vortons. // Each cell contains a dynamic array of integers // whose values are offsets into mVortons. UniformGrid <IntList> ugVortRef = new UniformGrid <IntList>(mInfluenceTree[0]); ugVortRef.Init(); int numVortons = mVortons.Count; for (int offset = 0 /* Start at 0th vorton */; offset < numVortons; ++offset) { // For each vorton... // Insert the vorton's offset into the spatial partition. ugVortRef[mVortons[offset].position].list.Add(offset); } // Phase 2: Exchange vorticity with nearest neighbors uint nx = ugVortRef.GetNumPoints(0); uint nxm1 = nx - 1; uint ny = ugVortRef.GetNumPoints(1); uint nym1 = ny - 1; uint nxy = nx * ny; uint nz = ugVortRef.GetNumPoints(2); uint nzm1 = nz - 1; uint[] idx = new uint[3]; for (idx[2] = 0; idx[2] < nzm1; ++idx[2]) { // For all points along z except the last... uint offsetZ0 = idx[2] * nxy; uint offsetZp = (idx[2] + 1) * nxy; for (idx[1] = 0; idx[1] < nym1; ++idx[1]) { // For all points along y except the last... uint offsetY0Z0 = idx[1] * nx + offsetZ0; uint offsetYpZ0 = (idx[1] + 1) * nx + offsetZ0; uint offsetY0Zp = idx[1] * nx + offsetZp; for (idx[0] = 0; idx[0] < nxm1; ++idx[0]) { // For all points along x except the last... uint offsetX0Y0Z0 = idx[0] + offsetY0Z0; for (int ivHere = 0; ivHere < ugVortRef[offsetX0Y0Z0].list.Count; ++ivHere) { // For each vorton in this gridcell... int rVortIdxHere = ugVortRef[offsetX0Y0Z0].list[ivHere]; Vorton rVortonHere = mVortons[rVortIdxHere]; Vector3 rVorticityHere = rVortonHere.vorticity; // Diffuse vorticity with other vortons in this same cell: for (int ivThere = ivHere + 1; ivThere < ugVortRef[offsetX0Y0Z0].list.Count; ++ivThere) { // For each OTHER vorton within this same cell... int rVortIdxThere = ugVortRef[offsetX0Y0Z0].list[ivThere]; Vorton rVortonThere = mVortons[rVortIdxThere]; Vector3 rVorticityThere = rVortonThere.vorticity; Vector3 vortDiff = rVorticityHere - rVorticityThere; Vector3 exchange = 2.0f * mViscosity * timeStep * vortDiff; // Amount of vorticity to exchange between particles. mVortons[rVortIdxHere].vorticity -= exchange; // Make "here" vorticity a little closer to "there". mVortons[rVortIdxThere].vorticity += exchange; // Make "there" vorticity a little closer to "here". } // Diffuse vorticity with vortons in adjacent cells: { uint offsetXpY0Z0 = idx[0] + 1 + offsetY0Z0; // offset of adjacent cell in +X direction for (int ivThere = 0; ivThere < ugVortRef[offsetXpY0Z0].list.Count; ++ivThere) { // For each vorton in the adjacent cell in +X direction... int rVortIdxThere = ugVortRef[offsetXpY0Z0].list[ivThere]; Vorton rVortonThere = mVortons[rVortIdxThere]; Vector3 rVorticityThere = rVortonThere.vorticity; Vector3 vortDiff = rVorticityHere - rVorticityThere; Vector3 exchange = mViscosity * timeStep * vortDiff; // Amount of vorticity to exchange between particles. mVortons[rVortIdxHere].vorticity -= exchange; // Make "here" vorticity a little closer to "there". mVortons[rVortIdxThere].vorticity += exchange; // Make "there" vorticity a little closer to "here". } } { uint offsetX0YpZ0 = idx[0] + offsetYpZ0; // offset of adjacent cell in +Y direction for (int ivThere = 0; ivThere < ugVortRef[offsetX0YpZ0].list.Count; ++ivThere) { // For each vorton in the adjacent cell in +Y direction... int rVortIdxThere = ugVortRef[offsetX0YpZ0].list[ivThere]; Vorton rVortonThere = mVortons[rVortIdxThere]; Vector3 rVorticityThere = rVortonThere.vorticity; Vector3 vortDiff = rVorticityHere - rVorticityThere; Vector3 exchange = mViscosity * timeStep * vortDiff; // Amount of vorticity to exchange between particles. mVortons[rVortIdxHere].vorticity -= exchange; // Make "here" vorticity a little closer to "there". mVortons[rVortIdxThere].vorticity += exchange; // Make "there" vorticity a little closer to "here". } } { uint offsetX0Y0Zp = idx[0] + offsetY0Zp; // offset of adjacent cell in +Z direction for (int ivThere = 0; ivThere < ugVortRef[offsetX0Y0Zp].list.Count; ++ivThere) { // For each vorton in the adjacent cell in +Z direction... int rVortIdxThere = ugVortRef[offsetX0Y0Zp].list[ivThere]; Vorton rVortonThere = mVortons[rVortIdxThere]; Vector3 rVorticityThere = rVortonThere.vorticity; Vector3 vortDiff = rVorticityHere - rVorticityThere; Vector3 exchange = mViscosity * timeStep * vortDiff; // Amount of vorticity to exchange between particles. mVortons[rVortIdxHere].vorticity -= exchange; // Make "here" vorticity a little closer to "there". mVortons[rVortIdxThere].vorticity += exchange; // Make "there" vorticity a little closer to "here". } } // Dissipate vorticity. See notes in header comment. mVortons[rVortIdxHere].vorticity -= mViscosity * timeStep * mVortons[rVortIdxHere].vorticity; // Reduce "here" vorticity. } } } } }