@@ -21,44 +21,74 @@ The computed solution for the default parameters looks like this:
2121module Example250_NSELidDrivenCavity
2222
2323using ExtendableFEM
24- using GridVisualize
24+ using Triangulate
2525using ExtendableGrids
26+ using SimplexGridFactory
2627using LinearAlgebra
28+ using GridVisualize
2729using Test # hide
2830
29- function kernel_nonlinear! (result, u_ops, qpinfo)
30- u, ∇u, p = view (u_ops, 1 : 2 ), view (u_ops, 3 : 6 ), view (u_ops, 7 )
31+
32+ function create_cone (h)
33+ builder = SimplexGridBuilder (; Generator = Triangulate)
34+
35+ # # points
36+ p1 = point! (builder, - 1 , 0 )
37+ p2 = point! (builder, 1 , 0 )
38+ p3 = point! (builder, 0 , - 2 )
39+
40+ # # top face
41+ facetregion! (builder, 1 )
42+ facet! (builder, p1, p2)
43+
44+ # # other faces
45+ facetregion! (builder, 2 )
46+ facet! (builder, p2, p3)
47+ facet! (builder, p3, p1)
48+
49+ cellregion! (builder, 1 )
50+ maxvolume! (builder, h)
51+ regionpoint! (builder, 0 , - 0.5 )
52+
53+ return simplexgrid (builder)
54+ end
55+
56+ const 𝕀 = [1 0 ; 0 1 ]
57+
58+ function NSE_kernel! (result, u_ops, qpinfo)
59+
60+ u = tensor_view (u_ops, 1 , TDVector (2 ))
61+ v = tensor_view (result, 1 , TDVector (2 ))
62+ ∇u = tensor_view (u_ops, 3 , TDMatrix (2 ))
63+ ∇v = tensor_view (result, 3 , TDMatrix (2 ))
64+ p = tensor_view (u_ops, 7 , TDScalar ())
65+ q = tensor_view (result, 7 , TDScalar ())
3166 μ = qpinfo. params[1 ]
32- result[1 ] = dot (u, view (∇u, 1 : 2 ))
33- result[2 ] = dot (u, view (∇u, 3 : 4 ))
34- result[3 ] = μ * ∇u[1 ] - p[1 ]
35- result[4 ] = μ * ∇u[2 ]
36- result[5 ] = μ * ∇u[3 ]
37- result[6 ] = μ * ∇u[4 ] - p[1 ]
38- result[7 ] = - (∇u[1 ] + ∇u[4 ])
67+
68+ tmul! (v, ∇u, u)
69+ ∇v .= μ .* ∇u .- p[1 ] .* 𝕀
70+ q[1 ] = - dot (∇u, 𝕀)
71+
3972 return nothing
4073end
4174
42- function boundarydata! (result, qpinfo)
75+
76+ # # boundary function
77+ function u_boundary! (result, qpinfo)
4378 result[1 ] = 1
4479 result[2 ] = 0
4580 return nothing
4681end
4782
48- function initialgrid_cone ()
49- xgrid = ExtendableGrid {Float64, Int32} ()
50- xgrid[Coordinates] = Array {Float64, 2} ([- 1 0 ; 0 - 2 ; 1 0 ]' )
51- xgrid[CellNodes] = Array {Int32, 2} ([1 2 3 ]' )
52- xgrid[CellGeometries] = VectorOfConstants {ElementGeometries, Int} (Triangle2D, 1 )
53- xgrid[CellRegions] = ones (Int32, 1 )
54- xgrid[BFaceRegions] = Array {Int32, 1} ([1 , 2 , 3 ])
55- xgrid[BFaceNodes] = Array {Int32, 2} ([1 2 ; 2 3 ; 3 1 ]' )
56- xgrid[BFaceGeometries] = VectorOfConstants {ElementGeometries, Int} (Edge1D, 3 )
57- xgrid[CoordinateSystem] = Cartesian2D
58- return xgrid
59- end
6083
61- function main (; μ_final = 0.0005 , order = 2 , nrefs = 5 , Plotter = nothing , kwargs... )
84+ function main (;
85+ μ_final = 0.0005 , # flow parameter
86+ order = 2 , # FE order of the flow field (pressure order is order-1)
87+ h = 1.0e-3 , # grid cell volume
88+ nrefs = 1 , # additional grid refinements
89+ Plotter = nothing ,
90+ kwargs...
91+ )
6292
6393 # # prepare parameter field
6494 extra_params = Array {Float64, 1} ([max (μ_final, 0.05 )])
@@ -70,17 +100,20 @@ function main(; μ_final = 0.0005, order = 2, nrefs = 5, Plotter = nothing, kwar
70100
71101 assign_unknown! (PD, u)
72102 assign_unknown! (PD, p)
73- assign_operator! (PD, NonlinearOperator (kernel_nonlinear!, [id (u), grad (u), id (p)]; params = extra_params, kwargs... ))
74- assign_operator! (PD, InterpolateBoundaryData (u, boundarydata!; regions = 3 ))
75- assign_operator! (PD, HomogeneousBoundaryData (u; regions = [1 , 2 ]))
76- assign_operator! (PD, FixDofs (p; dofs = [1 ]))
103+ assign_operator! (PD, NonlinearOperator (NSE_kernel!, [id (u), grad (u), id (p)]; params = extra_params, kwargs... ))
77104
105+ # # boundary data
106+ assign_operator! (PD, InterpolateBoundaryData (u, u_boundary!; regions = 1 ))
107+ assign_operator! (PD, HomogeneousBoundaryData (u; regions = [2 ]))
78108
79109 # # grid
80- xgrid = uniform_refine (initialgrid_cone ( ), nrefs)
110+ xgrid = uniform_refine (create_cone (h ), nrefs)
81111
82112 # # prepare FESpace
83- FES = [FESpace {H1Pk{2, 2, order}} (xgrid), FESpace {H1Pk{1, 2, order - 1}} (xgrid)]
113+ FES = [
114+ FESpace {H1Pk{2, 2, order}} (xgrid),
115+ FESpace {H1Pk{1, 2, order - 1}} (xgrid),
116+ ]
84117
85118 # # prepare plots
86119 plt = GridVisualizer (; Plotter = Plotter, layout = (1 , 2 ), clear = true , size = (1600 , 800 ))
@@ -93,7 +126,16 @@ function main(; μ_final = 0.0005, order = 2, nrefs = 5, Plotter = nothing, kwar
93126 while (true )
94127 step += 1
95128 @info " Step $step : solving for μ=$(extra_params[1 ]) "
96- sol, SC = ExtendableFEM. solve (PD, FES, SC; return_config = true , target_residual = 1.0e-10 , maxiterations = 20 , kwargs... )
129+ sol, SC = ExtendableFEM. solve (
130+ PD,
131+ FES,
132+ SC;
133+ return_config = true ,
134+ target_residual = 1.0e-10 ,
135+ maxiterations = 20 ,
136+ kwargs...
137+ )
138+
97139 if step == 1
98140 initialize! (PE, sol)
99141 end
@@ -116,9 +158,10 @@ function main(; μ_final = 0.0005, order = 2, nrefs = 5, Plotter = nothing, kwar
116158end
117159
118160generateplots = ExtendableFEM. default_generateplots (Example250_NSELidDrivenCavity, " example250.png" ) # hide
119- function runtests () # hide
120- sol, plt = main (; nrefs = 3 , μ_final = 0.005 ) # hide
121- @test sum (view (sol[1 ])) ≈ 9.501630403050289 # hide
122- return nothing # hide
123- end # hide
161+ function runtests () # hide
162+ sol, plt = main (; μ_final = 0.005 ) # hide
163+ sum (view (sol[1 ])) ≈ 237.24628017878518 # hide
164+ return nothing # hide
165+ end # hide
166+
124167end # module
0 commit comments