Nektar++
ImplicitExtrapolate.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: ImplicitExtrapolate.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: Abstract base class for ImplicitExtrapolate.
32//
33///////////////////////////////////////////////////////////////////////////////
34
36
37namespace Nektar
38{
39/**
40 * Registers the class with the Factory.
41 */
44 "Implicit", ImplicitExtrapolate::create, "Implicit");
45
48 eImplicit);
53 const SolverUtils::AdvectionSharedPtr advObject)
54 : WeakPressureExtrapolate(pSession, pFields, pPressure, pVel, advObject)
55{
56}
57
59{
60}
61
62/**
63 * Analogous to WeakPressureExtrapolate Implementation
64 * However, does not use Extrapolation of CurlCurl-term
65 * on boundary
66 */
68 const Array<OneD, const Array<OneD, NekDouble>> &fields,
69 const Array<OneD, const Array<OneD, NekDouble>> &N, NekDouble kinvis)
70{
72 if (m_HBCnumber > 0)
73 {
74 // Calculate just viscous BCs at current level and put in
75 // m_pressureHBCs[nlevels-1]
76 CalcNeumannPressureBCs(fields, N, kinvis);
77
78 // \int_bnd q n.u^{n} ds update current normal of field
79 // add m_pressureHBCs to gamma_0/Dt * m_acceleration[0]
80 AddVelBC();
81
82 // Copy m_pressureHBCs to m_PbndExp
84 }
85
86 // Evaluate High order outflow conditions
87 CalcOutflowBCs(fields, kinvis);
88}
89} // namespace Nektar
void CopyPressureHBCsToPbndExp(void)
void CalcNeumannPressureBCs(const Array< OneD, const Array< OneD, NekDouble > > &fields, const Array< OneD, const Array< OneD, NekDouble > > &N, NekDouble kinvis)
Definition: Extrapolate.h:172
int m_pressureCalls
number of times the high-order pressure BCs have been called
Definition: Extrapolate.h:232
void CalcOutflowBCs(const Array< OneD, const Array< OneD, NekDouble > > &fields, NekDouble kinvis)
static std::string solverTypeLookupId
ImplicitExtrapolate(const LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields, MultiRegions::ExpListSharedPtr pPressure, const Array< OneD, int > pVel, const SolverUtils::AdvectionSharedPtr advObject)
static std::string className
Name of class.
static ExtrapolateSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, MultiRegions::ExpListSharedPtr &pPressure, const Array< OneD, int > &pVel, const SolverUtils::AdvectionSharedPtr &advObject)
Creates an instance of this class.
void v_EvaluatePressureBCs(const Array< OneD, const Array< OneD, NekDouble > > &fields, const Array< OneD, const Array< OneD, NekDouble > > &N, NekDouble kinvis) override
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:197
static std::string RegisterEnumValue(std::string pEnum, std::string pString, int pEnumValue)
Registers an enumeration value.
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
std::shared_ptr< Advection > AdvectionSharedPtr
A shared pointer to an Advection object.
Definition: Advection.h:54
ExtrapolateFactory & GetExtrapolateFactory()
Definition: Extrapolate.cpp:48
double NekDouble