|
| 1 | +#! /usr/bin/env python3 |
| 2 | + |
| 3 | +# |
| 4 | +# Incompressible NSE solved within a channel geometry with parabolic inflow profile and an obstacle attached to the bottom towards the middle of the domain. The fluid field is initialized with a Stokes solution. The resulting velocity field is written to preCICE on the complete volume. |
| 5 | +# |
| 6 | + |
| 7 | +from nutils import function, mesh, cli, solver, export |
| 8 | +import treelog as log |
| 9 | +import numpy as np |
| 10 | +import precice |
| 11 | +from mpi4py import MPI |
| 12 | + |
| 13 | + |
| 14 | +def main(): |
| 15 | + |
| 16 | + print("Running Nutils") |
| 17 | + |
| 18 | + # define the Nutils mesh |
| 19 | + nx = 48 |
| 20 | + ny = 16 |
| 21 | + step_start = nx // 3 |
| 22 | + step_end = nx // 2 |
| 23 | + step_hight = ny // 2 |
| 24 | + |
| 25 | + grid = np.linspace(0, 6, nx + 1), np.linspace(0, 2, ny + 1) |
| 26 | + domain, geom = mesh.rectilinear(grid) |
| 27 | + domain = domain.withboundary(inflow="left", outflow="right", wall="top,bottom") - domain[ |
| 28 | + step_start:step_end, :step_hight |
| 29 | + ].withboundary(wall="left,top,right") |
| 30 | + |
| 31 | + # cloud of Gauss points |
| 32 | + gauss = domain.sample("gauss", degree=4) |
| 33 | + |
| 34 | + # Nutils namespace |
| 35 | + ns = function.Namespace() |
| 36 | + ns.x = geom |
| 37 | + |
| 38 | + ns.ubasis = domain.basis("std", degree=2).vector(2) |
| 39 | + ns.pbasis = domain.basis("std", degree=1) |
| 40 | + ns.u_i = "ubasis_ni ?u_n" # solution |
| 41 | + ns.p = "pbasis_n ?p_n" # solution |
| 42 | + ns.dudt_i = "ubasis_ni (?u_n - ?u0_n) / ?dt" # time derivative |
| 43 | + ns.μ = 0.5 # viscosity |
| 44 | + ns.σ_ij = "μ (u_i,j + u_j,i) - p δ_ij" |
| 45 | + ns.uin = "10 x_1 (2 - x_1)" # inflow profile |
| 46 | + |
| 47 | + # define the weak form, Stokes problem |
| 48 | + ures = gauss.integral("ubasis_ni,j σ_ij d:x" @ ns) |
| 49 | + pres = gauss.integral("pbasis_n u_k,k d:x" @ ns) |
| 50 | + |
| 51 | + # define Dirichlet boundary condition |
| 52 | + sqr = domain.boundary["inflow"].integral("(u_0 - uin)^2 d:x" @ ns, degree=2) |
| 53 | + sqr += domain.boundary["inflow,outflow"].integral("u_1^2 d:x" @ ns, degree=2) |
| 54 | + sqr += domain.boundary["wall"].integral("u_k u_k d:x" @ ns, degree=2) |
| 55 | + cons = solver.optimize(["u"], sqr, droptol=1e-15) |
| 56 | + |
| 57 | + # preCICE setup |
| 58 | + participant = precice.Participant("Fluid", "../precice-config.xml", 0, 1) |
| 59 | + |
| 60 | + # define coupling mesh |
| 61 | + mesh_name = "Fluid-Mesh" |
| 62 | + vertices = gauss.eval(ns.x) |
| 63 | + vertex_ids = participant.set_mesh_vertices(mesh_name, vertices) |
| 64 | + |
| 65 | + # coupling data |
| 66 | + data_name = "Velocity" |
| 67 | + |
| 68 | + participant.initialize() |
| 69 | + |
| 70 | + timestep = 0 |
| 71 | + solver_dt = 0.005 |
| 72 | + precice_dt = participant.get_max_time_step_size() |
| 73 | + dt = min(precice_dt, solver_dt) |
| 74 | + |
| 75 | + state = solver.solve_linear(("u", "p"), (ures, pres), constrain=cons) # initial condition |
| 76 | + |
| 77 | + # add convective term and time derivative for Navier-Stokes |
| 78 | + ures += gauss.integral("ubasis_ni (dudt_i + μ (u_i u_j)_,j) d:x" @ ns) |
| 79 | + |
| 80 | + while participant.is_coupling_ongoing(): |
| 81 | + |
| 82 | + if timestep % 1 == 0: # visualize |
| 83 | + bezier = domain.sample("bezier", 2) |
| 84 | + x, u, p = bezier.eval(["x_i", "u_i", "p"] @ ns, **state) |
| 85 | + with log.add(log.DataLog()): |
| 86 | + export.vtk("Fluid_" + str(timestep), bezier.tri, x, u=u, p=p) |
| 87 | + |
| 88 | + precice_dt = participant.get_max_time_step_size() |
| 89 | + |
| 90 | + # potentially adjust non-matching timestep sizes |
| 91 | + dt = min(solver_dt, precice_dt) |
| 92 | + |
| 93 | + # solve Nutils timestep |
| 94 | + state["u0"] = state["u"] |
| 95 | + state["dt"] = dt |
| 96 | + state = solver.newton(("u", "p"), (ures, pres), constrain=cons, arguments=state).solve(1e-10) |
| 97 | + |
| 98 | + velocity_values = gauss.eval(ns.u, **state) |
| 99 | + participant.write_data(mesh_name, data_name, vertex_ids, velocity_values) |
| 100 | + |
| 101 | + # do the coupling |
| 102 | + participant.advance(dt) |
| 103 | + |
| 104 | + # advance variables |
| 105 | + timestep += 1 |
| 106 | + |
| 107 | + participant.finalize() |
| 108 | + |
| 109 | + |
| 110 | +if __name__ == "__main__": |
| 111 | + cli.run(main) |
0 commit comments