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