|
| 1 | +#= |
| 2 | +
|
| 3 | +# 251 : Navier--Stokes Lid-driven cavity with iterative solver |
| 4 | +([source code](@__SOURCE_URL__)) |
| 5 | +
|
| 6 | +This example is a simpler version of Example250 and uses an iterative solver GMRES with ILU preconditioner. |
| 7 | +
|
| 8 | + |
| 9 | +=# |
| 10 | + |
| 11 | +module Example251_NSEWithIterativeSolver |
| 12 | + |
| 13 | +using ExtendableFEM |
| 14 | +#using ExtendableFEMBase |
| 15 | +using Triangulate |
| 16 | +using ExtendableGrids |
| 17 | +using SimplexGridFactory |
| 18 | +using LinearAlgebra |
| 19 | +using GridVisualize |
| 20 | +using LinearSolve |
| 21 | +using IncompleteLU |
| 22 | +using Test #hide |
| 23 | + |
| 24 | + |
| 25 | +function create_cone(h = 1.0e-2) |
| 26 | + builder = SimplexGridBuilder(; Generator = Triangulate) |
| 27 | + |
| 28 | + ## points |
| 29 | + p1 = point!(builder, -1, 0) |
| 30 | + p2 = point!(builder, 1, 0) |
| 31 | + p3 = point!(builder, 0, -2) |
| 32 | + |
| 33 | + ## top face |
| 34 | + facetregion!(builder, 1) |
| 35 | + facet!(builder, p1, p2) |
| 36 | + |
| 37 | + ## other faces |
| 38 | + facetregion!(builder, 2) |
| 39 | + facet!(builder, p2, p3) |
| 40 | + facet!(builder, p3, p1) |
| 41 | + |
| 42 | + cellregion!(builder, 1) |
| 43 | + maxvolume!(builder, h) |
| 44 | + regionpoint!(builder, 0, -0.5) |
| 45 | + |
| 46 | + return simplexgrid(builder) |
| 47 | +end |
| 48 | + |
| 49 | +const 𝕀 = [1 0; 0 1] |
| 50 | + |
| 51 | +function NSE_kernel!(result, u_ops, qpinfo) |
| 52 | + |
| 53 | + u = tensor_view(u_ops, 1, TDVector(2)) |
| 54 | + v = tensor_view(result, 1, TDVector(2)) |
| 55 | + ∇u = tensor_view(u_ops, 3, TDMatrix(2)) |
| 56 | + ∇v = tensor_view(result, 3, TDMatrix(2)) |
| 57 | + p = tensor_view(u_ops, 7, TDScalar()) |
| 58 | + q = tensor_view(result, 7, TDScalar()) |
| 59 | + μ = qpinfo.params[1] |
| 60 | + |
| 61 | + tmul!(v, ∇u, u) |
| 62 | + ∇v .= μ .* ∇u .- p[1] .* 𝕀 |
| 63 | + q[1] = -dot(∇u, 𝕀) |
| 64 | + |
| 65 | + return nothing |
| 66 | +end |
| 67 | + |
| 68 | +## boundary function |
| 69 | +function u_boundary!(result, qpinfo) |
| 70 | + result[1] = 1 |
| 71 | + result[2] = 0 |
| 72 | + return nothing |
| 73 | +end |
| 74 | + |
| 75 | +function main(; |
| 76 | + params = [0.05], # viscosity |
| 77 | + order = 2, # FE order of the flow field (pressure order is order-1) |
| 78 | + h = 1.0e-3, # grid cell volume |
| 79 | + Plotter = nothing, |
| 80 | + ) |
| 81 | + |
| 82 | + ## create the grid |
| 83 | + conegrid = create_cone(h) |
| 84 | + |
| 85 | + ## problem description |
| 86 | + PD = ProblemDescription("NSE") |
| 87 | + u = Unknown("u"; name = "velocity") |
| 88 | + p = Unknown("p"; name = "pressure") |
| 89 | + assign_unknown!(PD, u) |
| 90 | + assign_unknown!(PD, p) |
| 91 | + |
| 92 | + assign_operator!(PD, NonlinearOperator(NSE_kernel!, [id(u), grad(u), id(p)]; params = params, parallel = true)) |
| 93 | + |
| 94 | + # boundary values |
| 95 | + assign_operator!(PD, HomogeneousBoundaryData(u; value = 1, mask = [1, 0], regions = [1])) |
| 96 | + assign_operator!(PD, HomogeneousBoundaryData(u; value = 0, mask = [0, 1], regions = [1])) |
| 97 | + assign_operator!(PD, HomogeneousBoundaryData(u; regions = [2])) |
| 98 | + |
| 99 | + ## Taylor--Hood FE space |
| 100 | + FES = [ |
| 101 | + FESpace{H1Pk{2, 2, order}}(conegrid), |
| 102 | + FESpace{H1Pk{1, 2, order - 1}}(conegrid), |
| 103 | + ] |
| 104 | + |
| 105 | + ## prepare an initial solution matching the boundary data |
| 106 | + feVector = FEVector(FES; tags = PD.unknowns) |
| 107 | + interpolate!(feVector[u], ON_BFACES, u_boundary!; regions = [1]) |
| 108 | + |
| 109 | + ## solve with initial solution |
| 110 | + sol = ExtendableFEM.solve( |
| 111 | + PD, |
| 112 | + FES; |
| 113 | + init = feVector, |
| 114 | + method_linear = KrylovJL_GMRES(verbose = 1), |
| 115 | + precon_linear = IncompleteLU.ilu, # preconditioner! |
| 116 | + maxiterations = 10, |
| 117 | + ) |
| 118 | + |
| 119 | + ## plot streamlines |
| 120 | + PE = PointEvaluator([id(1)]) |
| 121 | + initialize!(PE, sol) |
| 122 | + plt = GridVisualizer(; Plotter, layout = (1, 1), clear = true, size = (800, 800)) |
| 123 | + GridVisualize.streamplot!(plt[1, 1], conegrid, eval_func_bary(PE), rasterpoints = 50, density = 2, title = "streamlines") |
| 124 | + |
| 125 | + return sol, plt |
| 126 | +end |
| 127 | + |
| 128 | +generateplots = ExtendableFEM.default_generateplots(Example251_NSEWithIterativeSolver, "example251.png") #hide |
| 129 | +function runtests() #hide |
| 130 | + sol, plt = main() #hide |
| 131 | + @test sum(view(sol[1])) ≈ 86.3125486317652 #hide |
| 132 | + return nothing #hide |
| 133 | +end #hide |
| 134 | + |
| 135 | +end # module |
0 commit comments