This notebook steps you through setting up fluid simulations and using TensorFlow's differentiation to optimize them.
Execute the cell below to install the ΦFlow Python package from GitHub.
# !pip install --quiet phiflow
from phi.flow import *
ΦFlow is vectorized but object-oriented, i.e. data are represented by Python objects that internally use tensors.
First, we create grids for the quantities we want to simulate. For this example, we require a velocity field and a smoke density field. We sample the smoke field at the cell centers and the velocity in staggered form.
smoke = CenteredGrid(0, extrapolation.BOUNDARY, x=32, y=40, bounds=Box(x=32, y=40)) # sampled at cell centers
velocity = StaggeredGrid(0, extrapolation.ZERO, x=32, y=40, bounds=Box(x=32, y=40)) # sampled in staggered form at face centers
Additionally, we want to add more smoke every time step.
We create the INFLOW
field from a circle (2D Sphere
) which defines where hot smoke is emitted.
Furthermore, we are interested in running the simulation for different inflow locations.
ΦFlow supports data-parallell execution via batch dimensions.
When a quantity has a batch dimension with size n, operations involving that quantity will be performed n times simultaneously and the result will also have that batch dimension. Here we add the batch dimension inflow_loc
.
For an overview of the dimension types, see the documentation or watch the introductory tutorial video.
INFLOW_LOCATION = tensor([(4, 5), (8, 5), (12, 5), (16, 5)], batch('inflow_loc'), channel(vector='x,y'))
INFLOW = 0.6 * CenteredGrid(Sphere(center=INFLOW_LOCATION, radius=3), extrapolation.BOUNDARY, x=32, y=40, bounds=Box(x=32, y=40))
The created grids are instances of the class Grid
.
Like tensors, grids also have the shape
attribute which lists all batch, spatial and channel dimensions.
Shapes in ΦFlow store not only the sizes of the dimensions but also their names and types.
print(f"Smoke: {smoke.shape}")
print(f"Velocity: {velocity.shape}")
print(f"Inflow: {INFLOW.shape}")
print(f"Inflow, spatial only: {INFLOW.shape.spatial}")
Smoke: (xˢ=32, yˢ=40) Velocity: (xˢ=32, yˢ=40, vectorᶜ=x,y) Inflow: (inflow_locᵇ=4, xˢ=32, yˢ=40) Inflow, spatial only: (xˢ=32, yˢ=40)
The grid values can be accessed using the values
property.
print(smoke.values)
print(velocity.values)
print(INFLOW.values)
(xˢ=32, yˢ=40) const 0.0 (~vectorᵈ=x,y, xˢ=~(x=31, y=32) int64, yˢ=~(x=40, y=39) int64) const 0.0 (inflow_locᵇ=4, xˢ=32, yˢ=40) 0.015 ± 0.094 (0e+00...6e-01)
Grids have many more properties which are documented here. Also note that the staggered grid has a non-uniform shape because the number of faces is not equal to the number of cells.
Next, let's do some physics!
Since the initial velocity is zero, we just add the inflow and the corresponding buoyancy force.
For the buoyancy force we use the factor (0, 0.5)
to specify strength and direction.
Finally, we project the velocity field to make it incompressible.
Note that the @
operator is a shorthand for resampling a field at different points. Since smoke
is sampled at cell centers and velocity
at face centers, this conversion is necessary.
smoke += INFLOW
buoyancy_force = smoke * (0, 0.5) @ velocity
velocity += buoyancy_force
velocity, _ = fluid.make_incompressible(velocity)
vis.plot(smoke)
/opt/hostedtoolcache/Python/3.12.7/x64/lib/python3.12/site-packages/phi/vis/_matplotlib/_matplotlib_plots.py:167: UserWarning: This figure includes Axes that are not compatible with tight_layout, so results might be incorrect. plt.tight_layout() # because subplot titles can be added after figure creation
Let's run a longer simulation! Now we add the transport or advection operations to the simulation. ΦFlow provides multiple algorithms for advection. Here we use semi-Lagrangian advection for the velocity and MacCormack advection for the smoke distribution.
trajectory = [smoke]
for i in range(20):
print(i, end=' ')
smoke = advect.mac_cormack(smoke, velocity, dt=1) + INFLOW
buoyancy_force = smoke * (0, 0.5) @ velocity
velocity = advect.semi_lagrangian(velocity, velocity, dt=1) + buoyancy_force
velocity, _ = fluid.make_incompressible(velocity)
trajectory.append(smoke)
trajectory = field.stack(trajectory, batch('time'))
vis.plot(trajectory, animate='time')
0 1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19
<Figure size 640x480 with 0 Axes>
The simulation we just computed was using pure NumPy so all operations were non-differentiable.
To enable differentiability, we need to use either PyTorch, TensorFlow or Jax.
This can be achieved by changing the import statement to phi.tf.flow
, phi.torch.flow
or phi.jax.flow
, respectively.
Tensors created after this import will be allocated using PyTorch / TensorFlow / Jax and operations on these will be executed with the corresponding backend.
These operations can make use of a GPU through CUDA if your configuration supports it.
# from phi.jax.flow import *
from phi.torch.flow import *
# from phi.tf.flow import *
We set up the simulation as before.
INFLOW_LOCATION = tensor([(4, 5), (8, 5), (12, 5), (16, 5)], batch('inflow_loc'), channel(vector='x,y'))
INFLOW = 0.6 * CenteredGrid(Sphere(center=INFLOW_LOCATION, radius=3), extrapolation.BOUNDARY, x=32, y=40, bounds=Box(x=32, y=40))
We can verify that tensors are now backed by TensorFlow / PyTorch / Jax.
type(INFLOW.values.native(INFLOW.shape))
torch.Tensor
Note that tensors created with NumPy will keep using NumPy/SciPy operations unless a TensorFlow tensor is also passed to the same operation.
Let's look at how to get gradients from our simulation.
Say we want to optimize the initial velocities so that all simulations arrive at a final state that is similar to the right simulation where the inflow is located at (16, 5)
.
To achieve this, we define the loss function as $L = | D(s - s_r) |^2$ where $s$ denotes the smoke density and the function $D$ diffuses the difference to smoothen the gradients.
def simulate(smoke: CenteredGrid, velocity: StaggeredGrid):
for _ in range(20):
smoke = advect.mac_cormack(smoke, velocity, dt=1) + INFLOW
buoyancy_force = smoke * (0, 0.5) @ velocity
velocity = advect.semi_lagrangian(velocity, velocity, dt=1) + buoyancy_force
velocity, _ = fluid.make_incompressible(velocity)
loss = math.sum(field.l2_loss(diffuse.explicit(smoke - field.stop_gradient(smoke.inflow_loc[-1]), 1, 1, 10)))
return loss, smoke, velocity
Now it is important that the initial velocity has the inflow_loc
dimension before we record the gradients.
initial_smoke = CenteredGrid(0, extrapolation.BOUNDARY, x=32, y=40, bounds=Box(x=32, y=40))
initial_velocity = StaggeredGrid(math.zeros(batch(inflow_loc=4)), extrapolation.ZERO, x=32, y=40, bounds=Box(x=32, y=40))
Finally, we use gradient_function()
to obtain the gradient with respect to the initial velocity. Since the velocity is the second argument to the simulate()
function, we pass wrt=[1]
.
sim_grad = field.gradient(simulate, wrt='velocity', get_output=False)
The argument get_output=False
specifies that we are not interested in the actual output of the function. By setting it to True
, we would also get the loss value and the final simulation state.
To evaluate the gradient, we simply call the gradient function with the same arguments as we would call the simulation.
velocity_grad = sim_grad(initial_smoke, initial_velocity)
vis.plot(velocity_grad)
/opt/hostedtoolcache/Python/3.12.7/x64/lib/python3.12/site-packages/phiml/backend/torch/_torch_backend.py:800: UserWarning: Sparse CSR tensor support is in beta state. If you miss a functionality in the sparse tensor support, please submit a feature request to https://github.com/pytorch/pytorch/issues. (Triggered internally at ../aten/src/ATen/SparseCsrTensorImpl.cpp:53.) return torch.sparse_csr_tensor(row_pointers, column_indices, values, shape, device=values.device)