Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
UnsteadyInviscidBurger.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File UnsteadyInviscidBurger.cpp
4 //
5 // For more information, please see: http://www.nektar.info
6 //
7 // The MIT License
8 //
9 // Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10 // Department of Aeronautics, Imperial College London (UK), and Scientific
11 // Computing and Imaging Institute, University of Utah (USA).
12 //
13 // License for the specific language governing rights and limitations under
14 // Permission is hereby granted, free of charge, to any person obtaining a
15 // copy of this software and associated documentation files (the "Software"),
16 // to deal in the Software without restriction, including without limitation
17 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
18 // and/or sell copies of the Software, and to permit persons to whom the
19 // Software is furnished to do so, subject to the following conditions:
20 //
21 // The above copyright notice and this permission notice shall be included
22 // in all copies or substantial portions of the Software.
23 //
24 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
25 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
26 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
27 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
28 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
29 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
30 // DEALINGS IN THE SOFTWARE.
31 //
32 // Description: Unsteady inviscid Burger solve routines
33 //
34 ///////////////////////////////////////////////////////////////////////////////
35 
37 
38 namespace Nektar
39 {
40  string UnsteadyInviscidBurger::className = GetEquationSystemFactory().RegisterCreatorFunction("UnsteadyInviscidBurger", UnsteadyInviscidBurger::create, "Inviscid Burger equation");
41 
44  : UnsteadySystem(pSession)
45  {
46  }
47 
48  /**
49  * @brief Initialisation object for the inviscid Burger equation.
50  */
52  {
53  // Call to the initialisation object of UnsteadySystem
55 
56  // Define the normal velocity fields
57  if (m_fields[0]->GetTrace())
58  {
59  m_traceVn = Array<OneD, NekDouble>(GetTraceNpoints());
60  }
61 
62  // Type of advection class to be used
63  switch(m_projectionType)
64  {
65  // Continuous field
67  {
68  string advName;
69  m_session->LoadSolverInfo("AdvectionType", advName, "NonConservative");
71  m_advection->SetFluxVector (&UnsteadyInviscidBurger::GetFluxVector, this);
72  break;
73  }
74  // Discontinuous field
76  {
77  string advName;
78  string riemName;
79 
80  m_session->LoadSolverInfo("AdvectionType", advName, "WeakDG");
82  m_advection->SetFluxVector (&UnsteadyInviscidBurger::GetFluxVector, this);
83 
84  m_session->LoadSolverInfo("UpwindType", riemName, "Upwind");
87 
88  m_advection->SetRiemannSolver(m_riemannSolver);
89  m_advection->InitObject (m_session, m_fields);
90  break;
91  }
92  default:
93  {
94  ASSERTL0(false, "Unsupported projection type.");
95  break;
96  }
97  }
98 
99  // If explicit it computes RHS and PROJECTION for the time integration
101  {
104  }
105  // Otherwise it gives an error because (no implicit integration)
106  else
107  {
108  ASSERTL0(false, "Implicit unsteady Advection not set up.");
109  }
110  }
111 
112  /**
113  * @brief Inviscid Burger equation destructor.
114  */
116  {
117  }
118 
119  /**
120  * @brief Get the normal velocity for the inviscid Burger equation.
121  */
123  {
124  // Counter variable
125  int i;
126 
127  // Number of trace (interface) points
128  int nTracePts = GetTraceNpoints();
129 
130  // Number of solution points
131  int nSolutionPts = GetNpoints();
132 
133  // Number of fields (variables of the problem)
134  int nVariables = m_fields.num_elements();
135 
136  // Auxiliary variables to compute the normal velocity
137  Array<OneD, NekDouble> Fwd (nTracePts);
138  Array<OneD, Array<OneD, NekDouble> > physfield (nVariables);
139 
140  // Reset the normal velocity
141  Vmath::Zero(nTracePts, m_traceVn, 1);
142 
143  // The TimeIntegration Class does not update the physical values of the
144  // solution. It is thus necessary to transform back the coefficient into
145  // the physical space and save them in physfield to compute the normal
146  // advection velocity properly. However it remains a critical point.
147  for(i = 0; i < nVariables; ++i)
148  {
149  physfield[i] = Array<OneD, NekDouble>(nSolutionPts);
150  m_fields[i]->BwdTrans_IterPerExp(m_fields[i]->GetCoeffs(),
151  physfield[i]);
152  }
153 
154  /// Extract the physical values at the trace space
155  m_fields[0]->ExtractTracePhys(physfield[0], Fwd);
156 
157  /// Compute the normal velocity
158  for (i = 0; i < m_spacedim; ++i)
159  {
160  Vmath::Vvtvp(nTracePts,
161  m_traceNormals[i], 1,
162  Fwd, 1,
163  m_traceVn, 1,
164  m_traceVn, 1);
165 
166  Vmath::Smul(nTracePts, 0.5, m_traceVn, 1, m_traceVn, 1);
167  }
168  return m_traceVn;
169  }
170 
171  /**
172  * @brief Compute the right-hand side for the inviscid Burger equation.
173  *
174  * @param inarray Given fields.
175  * @param outarray Calculated solution.
176  * @param time Time.
177  */
179  const Array<OneD, const Array<OneD, NekDouble> >&inarray,
180  Array<OneD, Array<OneD, NekDouble> >&outarray,
181  const NekDouble time)
182  {
183  // Counter variable
184  int i;
185 
186  // Number of fields (variables of the problem)
187  int nVariables = inarray.num_elements();
188 
189  // Number of solution points
190  int nSolutionPts = GetNpoints();
191 
192  // !Useless variable for WeakDG and FR!
193  Array<OneD, Array<OneD, NekDouble> > advVel;
194 
195  // RHS computation using the new advection base class
196  m_advection->Advect(nVariables, m_fields, advVel, inarray, outarray);
197 
198  // Negate the RHS
199  for (i = 0; i < nVariables; ++i)
200  {
201  Vmath::Neg(nSolutionPts, outarray[i], 1);
202  }
203  }
204 
205  /**
206  * @brief Compute the projection for the inviscid Burger equation.
207  *
208  * @param inarray Given fields.
209  * @param outarray Calculated solution.
210  * @param time Time.
211  */
213  const Array<OneD, const Array<OneD, NekDouble> >&inarray,
214  Array<OneD, Array<OneD, NekDouble> >&outarray,
215  const NekDouble time)
216  {
217  // Counter variable
218  int i;
219 
220  // Number of variables of the problem
221  int nVariables = inarray.num_elements();
222 
223  // Set the boundary conditions
224  SetBoundaryConditions(time);
225 
226  // Switch on the projection type (Discontinuous or Continuous)
227  switch(m_projectionType)
228  {
229  // Discontinuous projection
231  {
232  // Number of quadrature points
233  int nQuadraturePts = GetNpoints();
234 
235  // Just copy over array
236  for(i = 0; i < nVariables; ++i)
237  {
238  Vmath::Vcopy(nQuadraturePts, inarray[i], 1, outarray[i], 1);
239  }
240  break;
241  }
242 
243  // Continuous projection
246  {
247  Array<OneD, NekDouble> coeffs(m_fields[0]->GetNcoeffs());
248 
249  for(i = 0; i < nVariables; ++i)
250  {
251  m_fields[i]->FwdTrans(inarray[i], coeffs);
252  m_fields[i]->BwdTrans_IterPerExp(coeffs, outarray[i]);
253  }
254  break;
255  }
256  default:
257  ASSERTL0(false, "Unknown projection scheme");
258  break;
259  }
260  }
261 
262  /**
263  * @brief Return the flux vector for the inviscid Burger equation.
264  *
265  * @param i Component of the flux vector to calculate.
266  * @param physfield Fields.
267  * @param flux Resulting flux.
268  */
270  const Array<OneD, Array<OneD, NekDouble> > &physfield,
271  Array<OneD, Array<OneD, Array<OneD, NekDouble> > > &flux)
272  {
273  const int nq = GetNpoints();
274 
275  for (int i = 0; i < flux.num_elements(); ++i)
276  {
277  for (int j = 0; j < flux[0].num_elements(); ++j)
278  {
279  Vmath::Vmul(nq, physfield[i], 1, physfield[i], 1,
280  flux[i][j], 1);
281  Vmath::Smul(nq, 0.5, flux[i][j], 1, flux[i][j], 1);
282  }
283  }
284  }
285 }
286