Nektar++
NonSmoothShockCapture.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: NonSmoothShockCapture.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: NonSmooth artificial diffusion for shock capture
32//
33///////////////////////////////////////////////////////////////////////////////
34
36
37using namespace std;
38
39namespace Nektar
40{
41
45 "NonSmooth artificial diffusion for shock capture.");
46
50 const int spacedim)
51 : ArtificialDiffusion(pSession, pFields, spacedim)
52{
53 m_session->LoadParameter("SensorOffset", m_offset, 1);
54
55 // init h/p
56 m_varConv->SetElmtMinHP(m_fields);
57}
58
59/**
60 *
61 */
63 const Array<OneD, Array<OneD, NekDouble>> &physfield,
65{
66 int nPts = m_fields[0]->GetTotPoints();
67 int nElements = m_fields[0]->GetExpSize();
68
69 // Determine the maximum wavespeed
70 Array<OneD, NekDouble> Lambda(nPts, 0.0);
71 Array<OneD, NekDouble> soundspeed(nPts, 0.0);
72 Array<OneD, NekDouble> absVelocity(nPts, 0.0);
73 m_varConv->GetSoundSpeed(physfield, soundspeed);
74 m_varConv->GetAbsoluteVelocity(physfield, absVelocity);
75 Vmath::Vadd(nPts, absVelocity, 1, soundspeed, 1, Lambda, 1);
76
77 // Compute sensor based on rho
78 Array<OneD, NekDouble> Sensor(nPts, 0.0);
79 m_varConv->GetSensor(m_fields[0], physfield, Sensor, mu, m_offset);
80
81 // Scale AV
83 for (int e = 0; e < nElements; e++)
84 {
85 int physOffset = m_fields[0]->GetPhys_Offset(e);
86 int nElmtPoints = m_fields[0]->GetExp(e)->GetTotPoints();
87
88 // Get max wavespeed per element
89 NekDouble LambdaElmt = 0.0;
90 LambdaElmt = Vmath::Vmax(nElmtPoints, tmp = Lambda + physOffset, 1);
91
92 // Scale viscosity by the maximum wave speed and h/p
93 LambdaElmt *= m_mu0 * m_varConv->GetElmtMinHP()[e];
94 Vmath::Smul(nElmtPoints, LambdaElmt, tmp = mu + physOffset, 1,
95 tmp = mu + physOffset, 1);
96 }
97}
98} // namespace Nektar
Encapsulates the artificial diffusion used in shock capture.
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array of fields.
LibUtilities::SessionReaderSharedPtr m_session
Session reader.
NekDouble m_mu0
Constant scaling.
VariableConverterSharedPtr m_varConv
Auxiliary object to convert variables.
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
static ArtificialDiffusionSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const int spacedim)
Creates an instance of this class.
static std::string className
Name of the class.
void v_GetArtificialViscosity(const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, NekDouble > &mu) override
NonSmoothShockCapture(const LibUtilities::SessionReaderSharedPtr &pSession, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const int spacedim)
std::shared_ptr< SessionReader > SessionReaderSharedPtr
ArtificialDiffusionFactory & GetArtificialDiffusionFactory()
Declaration of the artificial diffusion factory singleton.
double NekDouble
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
T Vmax(int n, const T *x, const int incx)
Return the maximum element in x – called vmax to avoid conflict with max.
Definition: Vmath.hpp:644
STL namespace.