Nektar++
CoupledSolver.h
Go to the documentation of this file.
1namespace Nektar
2{
3/**
4* \page pageCoupledSolver Direct solver (coupled approach)
5
6 Here, we will explain how to use the incompressible Navier-Stokes solver of
7Nektar++ using the direct solver. Its source code is located at
8<code>Nektar++/solvers/IncNavierStokesSolver/EquationSystems</code>. We will
9start by showing the formulation and demonstrating how to run a first example
10with the incompressible Navier-Stokes solver and briefly explain the options of
11the solver specified in the input file. After these introductory explanations,
12we will present some numerical results in order to demonstrate the capabilities
13and the accuracy of the solver.
14 - @ref sectionFormulation_cuopled
15 - @ref sectionRunning1example_cuopled
16 - @ref sectionInputOptions_cuopled
17 - @ref sectionResult_cuopled
18 - @ref sectionRef_cuopled
19
20 \section sectionFormulation_cuopled Formulation
21
22 We consider the weak form of the Stokes problem for the velocity field
23\f$\boldsymbol{v}=[u,v]^{T}\f$ and the pressure field \f$p\f$: \f[
24 (\nabla \phi,\nu \nabla \boldsymbol{v}) -
25(\nabla\cdot\phi,p)=(\phi,\boldsymbol{f}) \f] \f[ (q,\nabla \cdot
26\boldsymbol{v}) = 0 \f]
27
28 where \f$\boldsymbol{v},\phi \in \boldsymbol{V}\f$, \f$p,q \in
29\boldsymbol{W}\f$ and \f$\boldsymbol{V},\boldsymbol{W}\f$ are appropriate spaces
30for the velocity and the pressure system to satisfy the inf-sup condition. Using
31a matrix notation the discrete system may be written as: \f{displaymath} \left[
32\begin{array}{ccc}
33 A & D_b^T & B\\
34 D_b & 0 & D_i^T\\
35 B^T & D_i & C
36 \end{array}\right]
37 \left[ \begin{array}{c}
38 \boldsymbol{v_b}\\
39 p\\
40 \boldsymbol{v_i}
41 \end{array}\right] =
42 \left[ \begin{array}{c}
43 \boldsymbol{f_b}\\
44 0\\
45 \boldsymbol{f_i}
46 \end{array}\right]
47 \f}
48
49 where the components of \f$A,B\f$ and \f$C\f$ are
50\f$(\nabla\phi_b,\nu\nabla\boldsymbol{v_b})\f$,
51\f$(\nabla\phi_b,\nu\nabla\boldsymbol{v_i})\f$ and
52 \f$(\nabla\phi_i,\nu\nabla\boldsymbol{v_i})\f$ and the components \f$D_b\f$ and
53\f$D_i\f$ are \f$(q,\nabla\boldsymbol{v_b})\f$ and
54\f$(q,\nabla\boldsymbol{v_i})\f$. The indices \f$b\f$ and \f$i\f$ refer to the
55degrees of freedom on the elemental boundary and interior, respectively. In
56constructing the system we have lumped the contributions form each component of
57the velocity field into matrices \f$A,B\f$ and \f$C\f$. However we note that for
58a Newtonian fluid the contribution from each field is decoupled. Since the
59inetrior degrees of freedom of the velocity field do not overlap, the matrix
60\f$C\f$ is block diagonal and to take advantage of this structure we can
61statically condense out the \f$C\f$ matrix to obtain the system: \f{displaymath}
62 \left[ \begin{array}{ccc}
63 A-BC^{-1}B^T & D_b^T-BC^{-1}D_i & 0\\
64 D_b-D_i^TC^{-1}B^T & -D_i^TC^{-1}D_i & 0\\
65 B^T & D_i & C
66 \end{array}\right]
67 \left[ \begin{array}{c}
68 \boldsymbol{v_b}\\
69 p\\
70 \boldsymbol{v_i}
71 \end{array}\right] =
72 \left[ \begin{array}{c}
73 \boldsymbol{f_b} - BC^{-1}\boldsymbol{f_i}\\
74 -D_i^TC^{-1}\boldsymbol{f_i}\\
75 \boldsymbol{f_i}
76 \end{array}\right]
77 \f}
78
79 To extend the aboce Stokes solver to an unsteady Navier-Stokes solver we first
80introduce the unsteady term, \f$\boldsymbol{v_t}\f$, into the Stokes problem.
81 This has the principal effect of modifying the weak Laplacian operator
82\f$(\nabla\phi,\nu\nabla\boldsymbol{v})\f$ into a weak Helmholts operator
83 \f$(\nabla\phi,\nu\nabla\boldsymbol{v})-\lambda(\phi,\boldsymbol{v})\f$ where
84\f$\lambda\f$ depends on the time integration scheme. The second modification
85requires the explicit discretisation of the non-linear terms in a similar manner
86to the splitting scheme and this term is then introduced as the forcing term
87\f$\boldsymbol{f}\f$.
88
89\section sectionRunning1example_cuopled Running a first example
90
91 The folder
92<code>Nektar++/regressionTests/Solvers/IncNavierStokesSolver/InputFiles</code>
93contains several <code>*.xml</code> files. These are input files for the
94Navier-Stokes solver specifying the geometry (i.e. the mesh and the spectal/hp
95expansion), the parameters and boundary conditions. Further details on the
96structure of the input file can be found in @ref pageXML.
97
98 Now, lets try to run one of them with the coupled solver.
99
100 - Copy the input file, <code>Test_ChanFlow_LinNS_m8.xml</code> to the directory
101where the solver is compiled, e.g.
102 <code>Nektar++/solvers/IncNavierStokesSolver/build/IncNavierStokesSolver</code>.
103 - Run the code by typing
104 @code
105 ./IncNavierStokesSolver Test_ChanFlow_LinNS_m8.xml
106 @endcode
107
108 The solution should now have been written to the file
109<code>Test_ChanFlow_LinNS_m8.fld</code>. This file is formatted in the Nektar++
110output format. To visualise the solution, we can convert the fld-file into
111TecPlot, Gmsh or Tecplot file formats using the Post-processing tools in
112<code>Nektar++/utilities/builds/PostProcessing/</code>. Here, we will
113demonstrate how to convert the <code>fld</code>-file into TecPlot-file format.
114
115 We convert the <code>fld</code>-file into Tecplot-file format by typing
116 @code
117 ../../../utilities/builds/PostProcessing/FldToTecplot
118Test_ChanFlow_LinNS_m8.xml Test_ChanFlow_LinNS_m8.fld
119 @endcode
120
121 It will create <code>Test_ChanFlow_LinNS_m8.dat</code> which can be loaed in
122TecPlot.
123
124 \image html CoupChanFlow.png "Channel Flow (u-velocity component)"
125
126\section sectionInputOptions_cuopled Input Options
127
128 A detailed descriptions of the input file for Nektar++ can be found in @ref
129pageXML. For what concern the incompressible Navier-Stokes solver a typical set
130is:
131 @code
132 <SOLVERINFO>
133 <I PROPERTY="SolverType" VALUE="CoupledLinearisedNS"/>
134 <I PROPERTY="EQTYPE" VALUE="UnsteadyNavierStokes"/>
135 <I PROPERTY="AdvectionForm" VALUE="Convective"/>
136 <I PROPERTY="Projection" VALUE="Galerkin"/>
137 <I PROPERTY="TimeIntegrationMethod" VALUE="IMEXOrder1"/>
138 </SOLVERINFO>
139 @endcode
140
141 @code
142 <PARAMETERS>
143 <P> TimeStep = 0.001 <P>
144 <P> NumSteps = 100 <P>
145 <P> IO_CheckSteps = 10 <P>
146 <P> IO_InfoSteps = 10 <P>
147 <P> Kinvis = 1.0 <P>
148 </PARAMETERS>
149 @endcode
150
151 Futher options can be specified in the input file, as for the Steady Oseen
152Flow:
153 @code
154 <SOLVERINFO>
155 <I PROPERTY="SolverType" VALUE="CoupledLinearisedNS"/>
156 <I PROPERTY="EQTYPE" VALUE="SteadyOseen"/>
157 <I PROPERTY="Projection" VALUE="Galerkin"/>
158 </SOLVERINFO>
159 @endcode
160
161 @code
162 <PARAMETERS>
163 <P> Kinvis = 0.025 </P>
164 <P> LAMBDA = 0.9637405441957 </P>
165 </PARAMETERS>
166 @endcode
167
168 @code
169 <FUNCTION NAME="AdvectionField">
170 <E VAR="u" VALUE="(1-exp(-LAMBDA*x)*cos(2*PI*y))" />
171 <E VAR="v" VALUE="(-LAMBDA/(2*PI))*exp(-LAMBDA*x)*sin(2*PI*y)" />
172 </FUNCTION>
173 @endcode
174
175\section sectionResult_cuopled Numerical Results
176
177 \image html Oseen.png "Steady Oseen Flow (u-velocity component)"
178
179
180\section sectionRef_cuopled References
181 [1] S.J. Sherwin and M. Ainsworth: <i>Unsteady Navier-Stokes Solvers Using
182Hybrid Spectral/hp Element Methods</i>, Conference Paper, 2000.<br /> [2] R.
183Stenberg and M. Suri: <i>Mixed hp finite element methods for problems in
184elasticity and Stokes flows</i>, Numer. Math, 72, 367-389, 1996.<br /> [3] M.
185Ainsworth and S.J. Sherwin: <i>Domain decomposition preconditioners for p and hp
186finite element approximation of the Stokes equations</i>, Comput. Methods Appl.
187Mech. Engrg., 175, 243-266, 1999. <br /> [4] P. La Tallec and A. Patra:
188<i>Non-overlapping domain decomposition methods for adaptive hp approximations
189of the Stokes problem with discontiunous pressure field</i>, Comput. Methods
190Appl. Mech. Engrg., 145, 361-379, 1997.<br /> [5] G. E. Karniadakis, M. Israeli,
191and S. A. Orszag: <i>High-order splitting methodsfor the incompressible
192Navier-Stokes equations</i>, J. Comput. Phys., 97, 414-443, 1991.<br />
193*/
194
195}