Nektar++
VelocityCorrectionSchemeImplicit.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File VelocityCorrectionSchemeImplicit.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// Permission is hereby granted, free of charge, to any person obtaining a
14// copy of this software and associated documentation files (the "Software"),
15// to deal in the Software without restriction, including without limitation
16// the rights to use, copy, modify, merge, publish, distribute, sublicense,
17// and/or sell copies of the Software, and to permit persons to whom the
18// Software is furnished to do so, subject to the following conditions:
19//
20// The above copyright notice and this permission notice shall be included
21// in all copies or substantial portions of the Software.
22//
23// THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24// OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25// FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26// THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27// LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28// FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29// DEALINGS IN THE SOFTWARE.
30//
31// Description: Implicit Velocity Correction Scheme for the Incompressible
32// Navier Stokes equations similar to Dong and Shen 2010
33///////////////////////////////////////////////////////////////////////////////
34
39
40#include <boost/algorithm/string.hpp>
41
42using namespace std;
43
44namespace Nektar
45{
46using namespace MultiRegions;
47
50 "VCSImplicit", VCSImplicit::create);
51
53 LibUtilities::SessionReader::RegisterEnumValue("SolverType", "VCSImplicit",
55
56/**
57 * Constructor. No special instructions here.
58 * v_DoInitialise sets the scheme up.
59 *
60 * \param
61 * \param
62 */
65 : UnsteadySystem(pSession, pGraph),
66 VelocityCorrectionScheme(pSession, pGraph)
67{
68}
69
70/**
71 * Destructor
72 */
74{
75}
76
77/**
78 *
79 */
81{
83 SolverUtils::AddSummaryItem(s, "Splitting Scheme",
84 "Implicit velocity correction");
85
86 if (m_extrapolation->GetSubStepName().size())
87 {
88 SolverUtils::AddSummaryItem(s, "Substepping",
89 m_extrapolation->GetSubStepName());
90 }
91
92 string dealias = m_homogen_dealiasing ? "Homogeneous1D" : "";
94 {
95 dealias += (dealias == "" ? "" : " + ") + string("spectral/hp");
96 }
97 if (dealias != "")
98 {
99 SolverUtils::AddSummaryItem(s, "Dealiasing", dealias);
100 }
101
102 string smoothing = m_useSpecVanVisc ? "spectral/hp" : "";
103 if (smoothing != "")
104 {
106 {
108 s, "Smoothing-SpecHP",
109 "SVV (" + smoothing + " Exp Kernel(cut-off = " +
110 boost::lexical_cast<string>(m_sVVCutoffRatio) +
111 ", diff coeff = " +
112 boost::lexical_cast<string>(m_sVVDiffCoeff) + "))");
113 }
114 else
115 {
117 {
119 s, "Smoothing-SpecHP",
120 "SVV (" + smoothing + " Power Kernel (Power ratio =" +
121 boost::lexical_cast<string>(m_sVVCutoffRatio) +
122 ", diff coeff = " +
123 boost::lexical_cast<string>(m_sVVDiffCoeff) +
124 "*Uh/p))");
125 }
126 else
127 {
129 s, "Smoothing-SpecHP",
130 "SVV (" + smoothing + " DG Kernel (diff coeff = " +
131 boost::lexical_cast<string>(m_sVVDiffCoeff) +
132 "*Uh/p))");
133 }
134 }
135 }
136
138 {
140 s, "Smoothing-Homo1D",
141 "SVV (Homogeneous1D - Exp Kernel(cut-off = " +
142 boost::lexical_cast<string>(m_sVVCutoffRatioHomo1D) +
143 ", diff coeff = " +
144 boost::lexical_cast<string>(m_sVVDiffCoeffHomo1D) + "))");
145 }
146
148 {
150 s, "GJP Stab. Impl. ",
151 m_session->GetSolverInfo("GJPStabilisation"));
152 SolverUtils::AddSummaryItem(s, "GJP Stab. JumpScale", m_GJPJumpScale);
153
154 if (boost::iequals(m_session->GetSolverInfo("GJPStabilisation"),
155 "Explicit"))
156 {
158 s, "GJP Normal Velocity",
159 m_session->GetSolverInfo("GJPNormalVelocity"));
160 }
161 }
162}
163
164/**
165 *
166 */
167void VCSImplicit::v_DoInitialise(bool dumpInitialConditions)
168{
169 VelocityCorrectionScheme::v_DoInitialise(dumpInitialConditions);
170
171 // Get integration order for extrapolation
172 m_intOrder = m_intScheme->GetOrder();
173
174 // Initialise Advection velocity
176 for (int i = 0; i < m_velocity.size(); i++)
177 {
179 }
180
181 // Initialise storage for explicit advection terms
183 for (int i = 0; i < m_nConvectiveFields; i++)
184 {
185 m_advection[i] =
187 }
188
189 // Check advection velocity: Extrapolated or Updated scheme
190 // Default to Extrapolated scheme
191 m_session->MatchSolverInfo("AdvectionVelocity", "Extrapolated",
192 m_advectionVelocity, true);
193
194 // Extrapolated velocity (for 2nd order scheme)
195 // Only initiliase for Extrapolated scheme
197 {
198 m_extVel =
200 for (int i = 0; i < m_velocity.size(); i++)
201 {
203 for (int j = 0; j < m_intOrder; j++)
204 {
205 m_extVel[i][j] =
207 }
208 }
209 }
210
211 // Check whether multiple or single GlobalLinSys are
212 // required and when to delete them
213 // Initialise to always reset GlobalLinSys
216
217 // Check skew-symmetric advection operator
218 // Default to Convective operator = false
219 m_session->MatchSolverInfo("AdvectionOperator", "SkewSymmetric",
221}
222
223/**
224 * Computes the forcing term for the pressure solve in
225 * VCSImplicit::v_SolvePressure(). It uses the weak pressure forcing similar to
226 * VCSWeakPressure (coefficient space).
227 *
228 * @param fields Holds the BDF formula of previous solutions \f$\sum_q
229 * \frac{\alpha_q}{\gamma} [u^{n-q}, v^{n-q}, \ldots]\f$.
230 * @param Forcing Array for the forcing term of the pressure solve. May
231 * contain temporary values on input. On output, Forcing = \f$ \int_\Omega
232 * \nabla \phi \cdot \left( \frac{1}{\Delta t} \sum_q \alpha_q \mathbf{u}^{n-q}
233 * + \mathbf{f}^{n+1} - N(\mathbf{u})^n \right) \f$.
234 *
235 */
237 const Array<OneD, const Array<OneD, NekDouble>> &fields,
239{
240 int ncoeffs = m_pressure->GetNcoeffs();
241 int phystot = m_fields[0]->GetTotPoints();
242 int nvel = m_velocity.size();
243
244 // Add 1/dt \sum_q \alpha_q u^{n-q} - N(u^n)
245 for (int i = 0; i < nvel; i++)
246 {
247 Vmath::Svtvp(phystot, 1.0 / aii_Dt, fields[i], 1, m_advection[i], 1,
248 Forcing[i], 1);
249 }
250
251 // Add forcing terms
252 for (auto &x : m_forcing)
253 {
255 }
256
257 // Do Innerproduct with derivative base
258 m_pressure->IProductWRTDerivBase(Forcing, Forcing[0]);
259
260 // Negate to reverse negation in HelmSolve
261 Vmath::Neg(ncoeffs, Forcing[0], 1);
262}
263
264/**
265 * Computes the forcing term and advection velocity for ADR solves.
266 *
267 * @param inarray Holds the BDF formula of previous solutions \f$\sum_q
268 \frac{\alpha_q}{\gamma} [u^{n-q}, v^{n-q}, \ldots]\f$.
269 * @param Forcing Shared forcing array m_F. Used in pressure and velocity
270 solve. On input, holds temporary values from pressure solve. On output: Forcing
271 = \f$ - \frac{1}{\nu} \left( \frac{1}{\Delta t} \sum_q \alpha_q
272 \mathbf{u}^{n-q} - \nabla p^{n+1} + \mathbf{f}^{n+1} \right) \f$.
273
274 * @param aii_Dt \f$\frac{\Delta t}{\gamma}\f$.
275 *
276 * Additionally, the advection velocity \f$ \tilde{\mathbf{u}}^{n+1} \f$ is
277 computed and handed to the velocity solve via the array m_AdvVel. It is
278 computed as \f$ \tilde{\mathbf{u}}^{n+1} = \sum_q \frac{\alpha_q}{\gamma}
279 \mathbf{u}^{n-q} - \frac{\Delta t}{\gamma} \left[ \nabla p^{n+1} -
280 N(\mathbf{u})^n - \nu \nabla \times \omega^n + \mathbf{f}^{n+1} \right] \f$.
281 *
282 */
284 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
286{
287 int phystot = m_fields[0]->GetTotPoints();
288 int nvel = m_velocity.size();
289
290 // Update (PhysSpace) pressure to n+1
291 m_pressure->BwdTrans(m_pressure->GetCoeffs(), m_pressure->UpdatePhys());
292
293 // Compute gradient \nabla p^{n+1}
294 if (nvel == 2)
295 {
296 m_pressure->PhysDeriv(m_pressure->GetPhys(), Forcing[m_velocity[0]],
297 Forcing[m_velocity[1]]);
298 }
299 else
300 {
301 m_pressure->PhysDeriv(m_pressure->GetPhys(), Forcing[m_velocity[0]],
303 }
304
305 // Negate pressure gradient: -1 * \nabla p
306 for (int i = 0; i < nvel; ++i)
307 {
308 Vmath::Neg(phystot, Forcing[i], 1);
309 }
310
311 // Zero convective fields
312 for (int i = nvel; i < m_nConvectiveFields; ++i)
313 {
314 Vmath::Zero(phystot, Forcing[i], 1);
315 }
316
317 // Add forcing terms: += f^{n+1}
318 for (auto &x : m_forcing)
319 {
321 }
322
323 // Build Advection Velocity
324 for (int i = 0; i < nvel; ++i)
325 {
326 // Extrapolated Advection Velocity = u^*n+1
328 {
329 Vmath::Smul(phystot, 1.0 / m_diffCoeff[i],
330 m_extVel[i][m_intOrder - 1], 1, m_AdvVel[i], 1);
331 }
332 // Updated Advection Velocity = u^*n+1 - (...)
333 else
334 {
335 // Get velocity fields [u^n, v^n, ..] for curl of vorticity
337 for (int i = 0; i < nvel; i++)
338 {
339 velocity[i] = Array<OneD, NekDouble>(m_fields[i]->GetPhys());
340 }
341 // Curl of vorticity \nabla \times \nabla \times [u^n, v^n, ..]
342 m_fields[0]->CurlCurl(velocity, velocity);
343
344 // \frac{\Delta t}{\gamma}(-\nabla p + f) + inarray
345 Vmath::Svtvp(phystot, aii_Dt, Forcing[i], 1, inarray[i], 1,
346 m_AdvVel[i], 1);
347 // += \frac{\Delta t}{\gamma} -N(u)
348 Vmath::Svtvp(phystot, aii_Dt, m_advection[i], 1, m_AdvVel[i], 1,
349 m_AdvVel[i], 1);
350 // += -\frac{\Delta t}{\gamma} CurlCurl(u)
351 Vmath::Svtvp(phystot, -aii_Dt * m_diffCoeff[i], velocity[i], 1,
352 m_AdvVel[i], 1, m_AdvVel[i], 1);
353 // *= 1/\nu
354 Vmath::Smul(phystot, 1.0 / m_diffCoeff[i], m_AdvVel[i], 1,
355 m_AdvVel[i], 1);
356 }
357 }
358
359 // Build Forcing term
360 for (int i = 0; i < nvel; ++i)
361 {
362 // \frac{\gamma}{\Delta t} * inarray - \nabla p + f
363 Vmath::Svtvp(phystot, 1.0 / aii_Dt, inarray[i], 1, Forcing[i], 1,
364 Forcing[i], 1);
365
366 // *= -1/kinvis
367 Vmath::Smul(phystot, -1.0 / m_diffCoeff[i], Forcing[i], 1, Forcing[i],
368 1);
369 }
370}
371
372/**
373 * Solve pressure system via a Poisson problem with ContField::HelmSolve().
374 * Uses coefficient space forcing and hence sets the
375 * argument PhysSpaceForcing=false in ContFieldHelmSolve().
376 *
377 * @param Forcing See output description in
378 * VCSImplicit::v_SetUpPressureForcing().
379 *
380 */
382{
384 // Setup coefficient for equation, Lambda = 0 ensures Poisson instead of
385 // Helmholtz problem
387
388 // Solve Pressure Poisson Equation (with Weak Forcing)
389 m_pressure->HelmSolve(Forcing, m_pressure->UpdateCoeffs(), factors,
392 false);
393
394 // Add presure to outflow bc if using convective like BCs
395 m_extrapolation->AddPressureToOutflowBCs(m_kinvis);
396}
397
398/**
399 * Solve velocity system via the
400 * ContField::LinearAdvectionDiffusionReactionSolve().
401 *
402 * @param Forcing Holds the forcing term for each velocity component. See output
403 * description in VCSImplicit::v_SetUpViscousForcing().
404 * @param inarray Unused in this routine. Still, holds \f$\sum_q
405 * \frac{\alpha_q}{\gamma} [u^{n-q}, v^{n-q}, \ldots] \f$.
406 * @param outarray Holds vector of previous solution \f$ [u^n, v^n, \ldots] \f$.
407 * Overwritten upon output with new solution [u^{n+1}, v^{n+1}, \ldots].
408 *
409 * Additionally, the advection velocity \f$ \tilde{\mathbf{u}}^{n+1} \f$ is
410 * handed to the velocity solve via the array m_AdvVel. It is computed in
411 * VCSImplicit::v_SetUpViscousForcing(). It is defined as \f$
412 * \tilde{\mathbf{u}}^{n+1} = \sum_q \frac{\alpha_q}{\gamma} \mathbf{u}^{n-q} -
413 * \frac{\Delta t}{\gamma} \left[ \nabla p^{n+1} - N(\mathbf{u})^n - \nu \nabla
414 * \times \omega^n + \mathbf{f}^{n+1} \right] \f$.
415 *
416 */
419 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
420 Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble aii_Dt)
421{
423 StdRegions::VarCoeffMap varcoeffs;
425
426 AppendSVVFactors(factors, varfactors);
427 ComputeGJPNormalVelocity(inarray, varcoeffs);
428
429 // Compute mass varcoeff for implicit Skew-symmetric operator
430 AddImplicitSkewSymAdvection(varcoeffs, aii_Dt);
431
432 // Set advection velocities
436 for (int i = 0; i < m_velocity.size(); i++)
437 {
438 varcoeffs[varcoefftypes[i]] = m_AdvVel[i];
439 }
440
441 // Solve Advection-Diffusion-Reaction system
442 for (int i = 0; i < m_nConvectiveFields; ++i)
443 {
444 // Add diffusion coefficient to GJP matrix operator (Implicit part)
446 {
448 }
449
450 // \lambda = - /frac{\gamma}{\Delta t \nu}
451 factors[StdRegions::eFactorLambda] = -1.0 / aii_Dt / m_diffCoeff[i];
452
453 auto gkey = m_fields[i]->LinearAdvectionDiffusionReactionSolve(
454 Forcing[i], m_fields[i]->UpdateCoeffs(), factors, varcoeffs,
455 varfactors);
456
457 // Nuke GlobalLinSys to avoid memory leak
458 // The condition nukes, if:
459 // - Remove after last velocity solve (w-velocity in 3D,
460 // assumes same matrix for each velocity)
461 // - Catches Null return from 3DH1D solve by checking matrix type
462 // is a LinearADR matrix (variant)
463 if (m_unsetGlobalLinSys[i] &&
464 (gkey.GetMatrixType() ==
466 gkey.GetMatrixType() ==
468 {
469 m_fields[i]->UnsetGlobalLinSys(gkey, true);
470 }
471
472 // Transform solution to PhysSpace
473 m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(), outarray[i]);
474 }
475}
476
477/**
478 * Explicit part of the method - HOPBCs based on VCSWeakPressure.
479 * For the implicit scheme, we do not extrapolate the advection and rotational
480 * diffusion terms, see also S. Dong and J. Shen (2010).
481 *
482 * @param inarray Holds the vector of previous solutions \f$ [u^n, v^n, \ldots]
483 * \f$.
484 * @param outarray In the implicit scheme, currently an empty array. TODO Don't
485 * waste memory here.
486 * @param time Holds the new time \f$ t^{n+1} = (n+1) \Delta t \f$.
487 *
488 */
490 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
491 Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble time)
492{
493 // Zero RHS derivatives, avoids undefined values
494 for (int i = 0; i < m_velocity.size(); ++i)
495 {
496 Vmath::Zero(outarray[i].size(), outarray[i], 1);
497 }
498
499 // Calculate High-Order pressure boundary conditions
501 timer.Start();
502 m_extrapolation->EvaluatePressureBCs(inarray, outarray, m_kinvis);
503 timer.Stop();
504 timer.AccumulateRegion("Pressure BCs");
505
506 // Extrapolate advection velocity with StifflyStableBeta coefficients
508 {
509 for (int i = 0; i < m_velocity.size(); ++i)
510 {
511 Vmath::Vcopy(inarray[i].size(), inarray[i], 1,
512 m_extVel[i][m_intOrder - 1], 1);
513 m_extrapolation->ExtrapolateArray(m_extVel[i]);
514 }
515 }
516
517 // Compute explicit advection terms
518 EvaluateAdvectionTerms(inarray, m_advection, time);
519}
520
522 NekDouble aii_Dt)
523{
525 {
526 return;
527 }
528
529 // Setup temporary storages
530 int nvel = m_velocity.size();
531 int phystot = m_fields[0]->GetTotPoints();
532 Array<OneD, NekDouble> divAdvVel(phystot, 0.0);
533 Array<OneD, NekDouble> tmpstore(phystot);
534 Array<OneD, NekDouble> ones(phystot, 1.0);
535
536 // Compute divergence of advection velocity: 1/\nu \tilde{u}
537 m_fields[0]->PhysDeriv(eX, m_AdvVel[0], divAdvVel);
538 for (int i = 1; i < nvel; ++i)
539 {
540 m_fields[i]->PhysDeriv(DirCartesianMap[i], m_AdvVel[i], tmpstore);
541 Vmath::Vadd(phystot, divAdvVel, 1, tmpstore, 1, divAdvVel, 1);
542 }
543 // Smul: 1/(2*lambda) * divAdvVel
544 NekDouble lambda = 1.0 / aii_Dt / m_diffCoeff[0];
545 Vmath::Smul(phystot, 1.0 / 2.0 / lambda, divAdvVel, 1, divAdvVel, 1);
546
547 // Vadd: 1 + 1/(2*lambda) * divAdvVel
548 Vmath::Vadd(phystot, ones, 1, divAdvVel, 1, divAdvVel, 1);
549
550 // Assign VarcoeffMass
551 varcoeffs[StdRegions::eVarCoeffMass] = divAdvVel;
552}
553
555{
556 // Loop over all convected Fields (excludes pressure)
557 for (int i = 0; i < m_nConvectiveFields; ++i)
558 {
559 // Always unset after last solve
560 if (i == m_nConvectiveFields - 1)
561 {
562 unsetGlobalLinSys[i] = 1;
563 continue;
564 }
565
566 auto contField =
567 std::dynamic_pointer_cast<MultiRegions::ContField>(m_fields[i]);
568 // Check against all solves after this one, whether same
569 // Global Linear system is required
570 for (int j = i + 1; j < m_nConvectiveFields; ++j)
571 {
572 // If same BCs (ie GlobalLinSys) in j-solves, save the GlobalLinSys
573 // for re-use
574 if (std::dynamic_pointer_cast<MultiRegions::ContField>(m_fields[j])
575 ->SameTypeOfBoundaryConditions(*contField))
576 {
577 unsetGlobalLinSys[i] = 0;
578 }
579 // If no same BCs (ie GlobalLinSys) in j-solves, UnsetGlobalLinSys()
580 // for i-solve
581 else
582 {
583 unsetGlobalLinSys[i] = 1;
584 }
585 }
586 }
587}
588
589} // namespace Nektar
MultiRegions::ExpListSharedPtr m_pressure
Pointer to field holding pressure field.
NekDouble m_kinvis
Kinematic viscosity.
ExtrapolateSharedPtr m_extrapolation
Array< OneD, int > m_velocity
int which identifies which components of m_fields contains the velocity (u,v,w);
int m_nConvectiveFields
Number of fields to be convected;.
std::vector< SolverUtils::ForcingSharedPtr > m_forcing
Forcing terms.
void EvaluateAdvectionTerms(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
static std::string RegisterEnumValue(std::string pEnum, std::string pString, int pEnumValue)
Registers an enumeration value.
void AccumulateRegion(std::string, int iolevel=0)
Accumulate elapsed time for a region.
Definition: Timer.cpp:70
NekDouble m_timestep
Time step size.
NekDouble m_time
Current time of simulation.
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
bool m_specHP_dealiasing
Flag to determine if dealisising is usde for the Spectral/hp element discretisation.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
enum HomogeneousType m_HomogeneousType
bool m_homogen_dealiasing
Flag to determine if dealiasing is used for homogeneous simulations.
SOLVER_UTILS_EXPORT int GetTotPoints()
Defines a forcing term to be explicitly applied.
Definition: Forcing.h:71
Base class for unsteady solvers.
LibUtilities::TimeIntegrationSchemeSharedPtr m_intScheme
Wrapper to the time integration scheme.
SOLVER_UTILS_EXPORT void v_GenerateSummary(SummaryList &s) override
Print a summary of time stepping parameters.
Array< OneD, NekInt > m_unsetGlobalLinSys
Array checking whether GlobalLinSys needs to be unset.
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_extVel
3D Array for extrapolated Advection Velocities [dir][time-levle][dof]
void v_SolveViscous(const Array< OneD, const Array< OneD, NekDouble > > &Forcing, const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble aii_Dt) override
void v_EvaluateAdvection_SetPressureBCs(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time) override
void v_SetUpPressureForcing(const Array< OneD, const Array< OneD, NekDouble > > &fields, Array< OneD, Array< OneD, NekDouble > > &Forcing, const NekDouble aii_Dt) override
void v_GenerateSummary(SolverUtils::SummaryList &s) override
Print a summary of time stepping parameters.
Array< OneD, Array< OneD, NekDouble > > m_AdvVel
void v_DoInitialise(bool dumpInitialConditions) override
Sets up initial conditions.
Array< OneD, Array< OneD, NekDouble > > m_advection
2D Array for Advection Velocities [dir][dof]
bool m_advectionVelocity
bool to identify implicit scheme ie advection velocity
void AddImplicitSkewSymAdvection(StdRegions::VarCoeffMap varcoeffs, NekDouble aii_Dt)
bool m_implicitSkewSymAdvection
bool to identify advection operator
VCSImplicit(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Constructor.
int m_intOrder
integer for advection velocity
static std::string className
Name of class.
void v_SetUpViscousForcing(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &Forcing, const NekDouble aii_Dt) override
void CheckUnsetGlobalLinSys(Array< OneD, NekInt > &unsetGlobalLinSys)
static SolverUtils::EquationSystemSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Creates an instance of this class.
void v_SolvePressure(const Array< OneD, NekDouble > &Forcing) override
Array< OneD, NekDouble > m_svvVarDiffCoeff
Array of coefficient if power kernel is used in SVV.
NekDouble m_sVVDiffCoeffHomo1D
Diffusion coefficient of SVV modes in homogeneous 1D Direction.
bool m_IsSVVPowerKernel
Identifier for Power Kernel otherwise DG kernel.
NekDouble m_sVVCutoffRatio
cutt off ratio from which to start decayhing modes
NekDouble m_sVVDiffCoeff
Diffusion coefficient of SVV modes.
bool m_useGJPStabilisation
bool to identify if GJP semi-implicit is active.
bool m_useSpecVanVisc
bool to identify if spectral vanishing viscosity is active.
Array< OneD, NekDouble > m_diffCoeff
Diffusion coefficients (will be kinvis for velocities)
void AppendSVVFactors(StdRegions::ConstFactorMap &factors, MultiRegions::VarFactorsMap &varFactorsMap)
void v_DoInitialise(bool dumpInitialConditions=true) override
Sets up initial conditions.
bool m_useHomo1DSpecVanVisc
bool to identify if spectral vanishing viscosity is active.
void ComputeGJPNormalVelocity(const Array< OneD, const Array< OneD, NekDouble > > &inarray, StdRegions::VarCoeffMap &varcoeffs)
std::shared_ptr< SessionReader > SessionReaderSharedPtr
MultiRegions::Direction const DirCartesianMap[]
Definition: ExpList.h:87
static VarFactorsMap NullVarFactorsMap
std::map< StdRegions::ConstFactorType, Array< OneD, NekDouble > > VarFactorsMap
std::vector< std::pair< std::string, std::string > > SummaryList
Definition: Misc.h:46
EquationSystemFactory & GetEquationSystemFactory()
void AddSummaryItem(SummaryList &l, const std::string &name, const std::string &value)
Adds a summary item to the summary info list.
Definition: Misc.cpp:47
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:174
@ eLinearAdvectionDiffusionReactionGJP
Definition: StdRegions.hpp:105
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:430
static VarCoeffMap NullVarCoeffMap
Definition: StdRegions.hpp:376
std::map< StdRegions::VarCoeffType, VarCoeffEntry > VarCoeffMap
Definition: StdRegions.hpp:375
StdRegions::ConstFactorMap factors
static Array< OneD, NekDouble > NullNekDouble1DArray
double NekDouble
void Svtvp(int n, const T alpha, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Svtvp (scalar times vector plus vector): z = alpha*x + y.
Definition: Vmath.hpp:396
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.hpp:292
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition: Vmath.hpp:180
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition: Vmath.hpp:100
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.hpp:273
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825
STL namespace.