Nektar++
Loading...
Searching...
No Matches
Static Public Member Functions | Static Public Attributes | Protected Member Functions | Protected Attributes | Private Member Functions | Private Attributes | Friends | List of all members
Nektar::LaplacePhi Class Reference

#include <LaplacePhi.h>

Inheritance diagram for Nektar::LaplacePhi:
[legend]

Static Public Member Functions

static EquationSystemSharedPtr create (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Creates an instance of this class.
 

Static Public Attributes

static std::string className
 Name of class.
 

Protected Member Functions

 LaplacePhi (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 
 ~LaplacePhi () override=default
 
void v_InitObject (bool DeclareFields=true) override
 Initialisation object for EquationSystem.
 
void v_GenerateSummary (SolverUtils::SummaryList &s) override
 Virtual function for generating summary information.
 
void v_DoSolve () override
 Virtual function for solve implementation.
 
void v_DoInitialise (bool dumpInitialConditions=true) override
 Virtual function for initialisation implementation.
 
- Protected Member Functions inherited from Nektar::SolverUtils::EquationSystem
SOLVER_UTILS_EXPORT EquationSystem (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Initialises EquationSystem class members.
 
virtual SOLVER_UTILS_EXPORT NekDouble v_LinfError (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray)
 Virtual function for the L_inf error computation between fields and a given exact solution.
 
virtual SOLVER_UTILS_EXPORT NekDouble v_L2Error (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray, bool Normalised=false)
 Virtual function for the L_2 error computation between fields and a given exact solution.
 
virtual SOLVER_UTILS_EXPORT NekDouble v_H1Error (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray, bool Normalised=false)
 Virtual function for the H_1 error computation between fields and a given exact solution.
 
virtual SOLVER_UTILS_EXPORT void v_TransCoeffToPhys ()
 Virtual function for transformation to physical space.
 
virtual SOLVER_UTILS_EXPORT void v_TransPhysToCoeff ()
 Virtual function for transformation to coefficient space.
 
virtual SOLVER_UTILS_EXPORT void v_SetInitialConditions (NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
 
virtual SOLVER_UTILS_EXPORT void v_EvaluateExactSolution (unsigned int field, Array< OneD, NekDouble > &outfield, const NekDouble time)
 
virtual SOLVER_UTILS_EXPORT void v_Output (void)
 
virtual SOLVER_UTILS_EXPORT MultiRegions::ExpListSharedPtr v_GetPressure (void)
 
virtual SOLVER_UTILS_EXPORT bool v_NegatedOp (void)
 Virtual function to identify if operator is negated in DoSolve.
 
virtual SOLVER_UTILS_EXPORT void v_ExtraFldOutput (std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
 

Protected Attributes

StdRegions::ConstFactorMap m_factors
 
- Protected Attributes inherited from Nektar::SolverUtils::EquationSystem
LibUtilities::CommSharedPtr m_comm
 Communicator.
 
bool m_verbose
 
LibUtilities::SessionReaderSharedPtr m_session
 The session reader.
 
std::map< std::string, SolverUtils::SessionFunctionSharedPtrm_sessionFunctions
 Map of known SessionFunctions.
 
LibUtilities::FieldIOSharedPtr m_fld
 Field input/output.
 
Array< OneD, MultiRegions::ExpListSharedPtrm_fields
 Array holding all dependent variables.
 
SpatialDomains::BoundaryConditionsSharedPtr m_boundaryConditions
 Pointer to boundary conditions object.
 
SpatialDomains::MeshGraphSharedPtr m_graph
 Pointer to graph defining mesh.
 
std::string m_sessionName
 Name of the session.
 
NekDouble m_time
 Current time of simulation.
 
int m_initialStep
 Number of the step where the simulation should begin.
 
NekDouble m_fintime
 Finish time of the simulation.
 
NekDouble m_timestep
 Time step size.
 
NekDouble m_lambda
 Lambda constant in real system if one required.
 
NekDouble m_checktime
 Time between checkpoints.
 
NekDouble m_lastCheckTime
 
NekDouble m_TimeIncrementFactor
 
int m_nchk
 Number of checkpoints written so far.
 
int m_steps
 Number of steps to take.
 
int m_checksteps
 Number of steps between checkpoints.
 
int m_infosteps
 Number of time steps between outputting status information.
 
int m_iterPIT = 0
 Number of parallel-in-time time iteration.
 
int m_windowPIT = 0
 Index of windows for parallel-in-time time iteration.
 
int m_spacedim
 Spatial dimension (>= expansion dim).
 
int m_expdim
 Expansion dimension.
 
bool m_singleMode
 Flag to determine if single homogeneous mode is used.
 
bool m_halfMode
 Flag to determine if half homogeneous mode is used.
 
bool m_multipleModes
 Flag to determine if use multiple homogenenous modes are used.
 
bool m_useFFT
 Flag to determine if FFT is used for homogeneous transform.
 
bool m_homogen_dealiasing
 Flag to determine if dealiasing is used for homogeneous simulations.
 
bool m_specHP_dealiasing
 Flag to determine if dealisising is usde for the Spectral/hp element discretisation.
 
enum MultiRegions::ProjectionType m_projectionType
 Type of projection; e.g continuous or discontinuous.
 
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
 Array holding trace normals for DG simulations in the forwards direction.
 
Array< OneD, bool > m_checkIfSystemSingular
 Flag to indicate if the fields should be checked for singularity.
 
LibUtilities::FieldMetaDataMap m_fieldMetaDataMap
 Map to identify relevant solver info to dump in output fields.
 
Array< OneD, NekDoublem_movingFrameData
 Moving reference frame status in the inertial frame X, Y, Z, Theta_x, Theta_y, Theta_z, U, V, W, Omega_x, Omega_y, Omega_z, A_x, A_y, A_z, DOmega_x, DOmega_y, DOmega_z, pivot_x, pivot_y, pivot_z.
 
std::vector< std::string > m_strFrameData
 variable name in m_movingFrameData
 
int m_NumQuadPointsError
 Number of Quadrature points used to work out the error.
 
enum HomogeneousType m_HomogeneousType
 
NekDouble m_LhomX
 physical length in X direction (if homogeneous)
 
NekDouble m_LhomY
 physical length in Y direction (if homogeneous)
 
NekDouble m_LhomZ
 physical length in Z direction (if homogeneous)
 
int m_npointsX
 number of points in X direction (if homogeneous)
 
int m_npointsY
 number of points in Y direction (if homogeneous)
 
int m_npointsZ
 number of points in Z direction (if homogeneous)
 
int m_HomoDirec
 number of homogenous directions
 

Private Member Functions

Array< OneD, bool > v_GetSystemSingularChecks () override
 
void setUserDefinedBC (int n)
 
void calculateAddedMass ()
 

Private Attributes

Array< OneD, Array< OneD, NekDouble > > m_boundValues
 
Array< OneD, NekDoublem_pivot
 

Friends

class MemoryManager< LaplacePhi >
 Class may only be instantiated through the MemoryManager.
 

Additional Inherited Members

- Public Member Functions inherited from Nektar::SolverUtils::EquationSystem
virtual SOLVER_UTILS_EXPORT ~EquationSystem ()
 Destructor.
 
SOLVER_UTILS_EXPORT void InitObject (bool DeclareField=true)
 Initialises the members of this object.
 
SOLVER_UTILS_EXPORT void DoInitialise (bool dumpInitialConditions=true)
 Perform any initialisation necessary before solving the problem.
 
SOLVER_UTILS_EXPORT void DoSolve ()
 Solve the problem.
 
SOLVER_UTILS_EXPORT void TransCoeffToPhys ()
 Transform from coefficient to physical space.
 
SOLVER_UTILS_EXPORT void TransPhysToCoeff ()
 Transform from physical to coefficient space.
 
SOLVER_UTILS_EXPORT void Output ()
 Perform output operations after solve.
 
SOLVER_UTILS_EXPORT std::string GetSessionName ()
 Get Session name.
 
template<class T >
std::shared_ptr< T > as ()
 
SOLVER_UTILS_EXPORT void ResetSessionName (std::string newname)
 Reset Session name.
 
SOLVER_UTILS_EXPORT LibUtilities::SessionReaderSharedPtr GetSession ()
 Get Session name.
 
SOLVER_UTILS_EXPORT MultiRegions::ExpListSharedPtr GetPressure ()
 Get pressure field if available.
 
SOLVER_UTILS_EXPORT void ExtraFldOutput (std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
 
SOLVER_UTILS_EXPORT void PrintSummary (std::ostream &out)
 Print a summary of parameters and solver characteristics.
 
SOLVER_UTILS_EXPORT void SetLambda (NekDouble lambda)
 Set parameter m_lambda.
 
SOLVER_UTILS_EXPORT SessionFunctionSharedPtr GetFunction (std::string name, const MultiRegions::ExpListSharedPtr &field=MultiRegions::NullExpListSharedPtr, bool cache=false)
 Get a SessionFunction by name.
 
SOLVER_UTILS_EXPORT void SetInitialConditions (NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
 Initialise the data in the dependent fields.
 
SOLVER_UTILS_EXPORT void EvaluateExactSolution (int field, Array< OneD, NekDouble > &outfield, const NekDouble time)
 Evaluates an exact solution.
 
SOLVER_UTILS_EXPORT NekDouble L2Error (unsigned int field, const Array< OneD, NekDouble > &exactsoln, bool Normalised=false)
 Compute the L2 error between fields and a given exact solution.
 
SOLVER_UTILS_EXPORT NekDouble L2Error (unsigned int field, bool Normalised=false)
 Compute the L2 error of the fields.
 
SOLVER_UTILS_EXPORT NekDouble LinfError (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray)
 Linf error computation.
 
SOLVER_UTILS_EXPORT NekDouble H1Error (unsigned int field, const Array< OneD, NekDouble > &exactsoln, bool Normalised=false)
 Compute the H1 error between fields and a given exact solution.
 
SOLVER_UTILS_EXPORT Array< OneD, NekDoubleErrorExtraPoints (unsigned int field)
 Compute error (L2 and L_inf) over an larger set of quadrature points return [L2 Linf].
 
SOLVER_UTILS_EXPORT void Checkpoint_Output (const int n)
 Write checkpoint file of m_fields.
 
SOLVER_UTILS_EXPORT void Checkpoint_Output (const int n, MultiRegions::ExpListSharedPtr &field, std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
 Write checkpoint file of custom data fields.
 
SOLVER_UTILS_EXPORT void Checkpoint_BaseFlow (const int n)
 Write base flow file of m_fields.
 
SOLVER_UTILS_EXPORT void WriteFld (const std::string &outname)
 Write field data to the given filename.
 
SOLVER_UTILS_EXPORT void WriteFld (const std::string &outname, MultiRegions::ExpListSharedPtr &field, std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
 Write input fields to the given filename.
 
SOLVER_UTILS_EXPORT void ImportFld (const std::string &infile, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields)
 Input field data from the given file.
 
SOLVER_UTILS_EXPORT void ImportFldToMultiDomains (const std::string &infile, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const int ndomains)
 Input field data from the given file to multiple domains.
 
SOLVER_UTILS_EXPORT void ImportFld (const std::string &infile, std::vector< std::string > &fieldStr, Array< OneD, Array< OneD, NekDouble > > &coeffs)
 Output a field. Input field data into array from the given file.
 
SOLVER_UTILS_EXPORT void ImportFld (const std::string &infile, MultiRegions::ExpListSharedPtr &pField, std::string &pFieldName)
 Output a field. Input field data into ExpList from the given file.
 
SOLVER_UTILS_EXPORT void SessionSummary (SummaryList &vSummary)
 Write out a session summary.
 
SOLVER_UTILS_EXPORT Array< OneD, MultiRegions::ExpListSharedPtr > & UpdateFields ()
 
SOLVER_UTILS_EXPORT LibUtilities::FieldMetaDataMapUpdateFieldMetaDataMap ()
 Get hold of FieldInfoMap so it can be updated.
 
SOLVER_UTILS_EXPORT NekDouble GetTime ()
 Return final time.
 
SOLVER_UTILS_EXPORT int GetNcoeffs ()
 
SOLVER_UTILS_EXPORT int GetNcoeffs (const int eid)
 
SOLVER_UTILS_EXPORT int GetNumExpModes ()
 
SOLVER_UTILS_EXPORT const Array< OneD, int > GetNumExpModesPerExp ()
 
SOLVER_UTILS_EXPORT int GetNvariables ()
 
SOLVER_UTILS_EXPORT const std::string GetVariable (unsigned int i)
 
SOLVER_UTILS_EXPORT int GetTraceTotPoints ()
 
SOLVER_UTILS_EXPORT int GetTraceNpoints ()
 
SOLVER_UTILS_EXPORT int GetExpSize ()
 
SOLVER_UTILS_EXPORT int GetPhys_Offset (int n)
 
SOLVER_UTILS_EXPORT int GetCoeff_Offset (int n)
 
SOLVER_UTILS_EXPORT int GetTotPoints ()
 
SOLVER_UTILS_EXPORT int GetTotPoints (int n)
 
SOLVER_UTILS_EXPORT int GetNpoints ()
 
SOLVER_UTILS_EXPORT int GetSteps ()
 
SOLVER_UTILS_EXPORT void SetSteps (const int steps)
 
SOLVER_UTILS_EXPORT NekDouble GetTimeStep ()
 
SOLVER_UTILS_EXPORT void CopyFromPhysField (const int i, Array< OneD, NekDouble > &output)
 
SOLVER_UTILS_EXPORT void CopyToPhysField (const int i, const Array< OneD, const NekDouble > &input)
 
SOLVER_UTILS_EXPORT Array< OneD, NekDouble > & UpdatePhysField (const int i)
 
SOLVER_UTILS_EXPORT void ZeroPhysFields ()
 
SOLVER_UTILS_EXPORT void FwdTransFields ()
 
SOLVER_UTILS_EXPORT void SetModifiedBasis (const bool modbasis)
 
SOLVER_UTILS_EXPORT int GetCheckpointNumber ()
 
SOLVER_UTILS_EXPORT void SetCheckpointNumber (int num)
 
SOLVER_UTILS_EXPORT int GetCheckpointSteps ()
 
SOLVER_UTILS_EXPORT void SetCheckpointSteps (int num)
 
SOLVER_UTILS_EXPORT int GetInfoSteps ()
 
SOLVER_UTILS_EXPORT void SetInfoSteps (int num)
 
SOLVER_UTILS_EXPORT void SetIterationNumberPIT (int num)
 
SOLVER_UTILS_EXPORT void SetWindowNumberPIT (int num)
 
SOLVER_UTILS_EXPORT Array< OneD, const Array< OneD, NekDouble > > GetTraceNormals ()
 
SOLVER_UTILS_EXPORT void SetTime (const NekDouble time)
 
SOLVER_UTILS_EXPORT void SetTimeStep (const NekDouble timestep)
 
SOLVER_UTILS_EXPORT void SetInitialStep (const int step)
 
SOLVER_UTILS_EXPORT void SetBoundaryConditions (NekDouble time)
 Evaluates the boundary conditions at the given time.
 
SOLVER_UTILS_EXPORT bool NegatedOp ()
 Identify if operator is negated in DoSolve.
 
- Protected Types inherited from Nektar::SolverUtils::EquationSystem
enum  HomogeneousType { eHomogeneous1D , eHomogeneous2D , eHomogeneous3D , eNotHomogeneous }
 Parameter for homogeneous expansions. More...
 
- Static Protected Attributes inherited from Nektar::SolverUtils::EquationSystem
static std::string equationSystemTypeLookupIds []
 
static std::string projectionTypeLookupIds []
 

Detailed Description

Definition at line 44 of file LaplacePhi.h.

Constructor & Destructor Documentation

◆ LaplacePhi()

Nektar::LaplacePhi::LaplacePhi ( const LibUtilities::SessionReaderSharedPtr pSession,
const SpatialDomains::MeshGraphSharedPtr pGraph 
)
protected

Definition at line 45 of file LaplacePhi.cpp.

47 : EquationSystem(pSession, pGraph), m_factors()
48{
51}
StdRegions::ConstFactorMap m_factors
Definition LaplacePhi.h:65
A base class for describing how to solve specific equations.

References Nektar::StdRegions::eFactorLambda, Nektar::StdRegions::eFactorTau, and m_factors.

◆ ~LaplacePhi()

Nektar::LaplacePhi::~LaplacePhi ( )
overrideprotecteddefault

Member Function Documentation

◆ calculateAddedMass()

void Nektar::LaplacePhi::calculateAddedMass ( )
private

Definition at line 203 of file LaplacePhi.cpp.

204{
205 Array<OneD, const SpatialDomains::BoundaryConditionShPtr> m_BndConds =
206 m_fields[0]->GetBndConditions();
207 Array<OneD, MultiRegions::ExpListSharedPtr> m_BndExp =
208 m_fields[0]->GetBndCondExpansions();
209 int nfld = m_fields.size();
210 Array<OneD, NekDouble> value(nfld * nfld,
211 std::numeric_limits<NekDouble>::lowest());
212 int offset = 0;
213 for (int i = 0; i < m_BndConds.size(); ++i)
214 {
215 if (boost::iequals(m_BndConds[i]->GetUserDefined(), "MOVE"))
216 {
217 MultiRegions::ExpListSharedPtr edgeExplist = m_BndExp[i];
218 Array<OneD, Array<OneD, NekDouble>> phi(nfld);
220 m_fields[0]->GetBndElmtExpansion(i, BndElmtExp, false);
221 Array<OneD, NekDouble> phiElm(BndElmtExp->GetTotPoints(), 0.);
222 for (int j = 0; j < nfld; ++j)
223 {
224 m_fields[j]->ExtractPhysToBndElmt(i, m_fields[j]->GetPhys(),
225 phiElm);
226 m_fields[j]->ExtractElmtToBndPhys(i, phiElm, phi[j]);
227 }
228 ////phi_0 n_0, phi_0 n_1, phi_1 n_0, phi_1 n_1
229 Array<OneD, NekDouble> mularray(edgeExplist->GetTotPoints(), 0.);
230 for (int j = 0; j < nfld; ++j)
231 {
232 for (int k = 0; k < m_boundValues.size(); ++k)
233 {
234 Vmath::Vmul(edgeExplist->GetTotPoints(), phi[j], 1,
235 m_boundValues[k] + offset, 1, mularray, 1);
236 value[k + j * nfld] = edgeExplist->Integral(mularray);
237 }
238 }
239 offset += edgeExplist->GetTotPoints();
240 }
241 }
242
243 m_session->GetComm()->AllReduce(value, LibUtilities::ReduceMax);
244 if (m_session->GetComm()->GetRank() == 0)
245 {
246 for (int j = 0; j < nfld; ++j)
247 {
248 for (int k = 0; k < nfld; ++k)
249 {
250 std::cout << "value[" << m_session->GetVariable(j)[3] << ", "
251 << m_session->GetVariable(k)[3]
252 << "] = " << std::scientific << std::setprecision(7)
253 << value[k + j * nfld] << std::endl;
254 }
255 }
256 }
257}
Array< OneD, Array< OneD, NekDouble > > m_boundValues
Definition LaplacePhi.h:81
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition Vmath.hpp:72

References m_boundValues, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_session, Nektar::LibUtilities::ReduceMax, and Vmath::Vmul().

Referenced by v_DoSolve().

◆ create()

static EquationSystemSharedPtr Nektar::LaplacePhi::create ( const LibUtilities::SessionReaderSharedPtr pSession,
const SpatialDomains::MeshGraphSharedPtr pGraph 
)
inlinestatic

Creates an instance of this class.

Definition at line 51 of file LaplacePhi.h.

54 {
57 p->InitObject();
58 return p;
59 }
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
std::shared_ptr< EquationSystem > EquationSystemSharedPtr
A shared pointer to an EquationSystem object.
std::vector< double > p(NPUPPER)

References Nektar::MemoryManager< DataType >::AllocateSharedPtr().

◆ setUserDefinedBC()

void Nektar::LaplacePhi::setUserDefinedBC ( int  n)
private

Definition at line 184 of file LaplacePhi.cpp.

185{
186 int m = m_session->GetVariable(n)[3] - '0';
187 Array<OneD, const SpatialDomains::BoundaryConditionShPtr> m_BndConds =
188 m_fields[n]->GetBndConditions();
189 Array<OneD, MultiRegions::ExpListSharedPtr> m_BndExp =
190 m_fields[n]->GetBndCondExpansions();
191 int offset = 0;
192 for (int i = 0; i < m_BndConds.size(); ++i)
193 {
194 if (boost::iequals(m_BndConds[i]->GetUserDefined(), "MOVE"))
195 {
196 m_BndExp[i]->IProductWRTBase(m_boundValues[m] + offset,
197 m_BndExp[i]->UpdateCoeffs());
198 offset += m_BndExp[i]->GetTotPoints();
199 }
200 }
201}

References m_boundValues, Nektar::SolverUtils::EquationSystem::m_fields, and Nektar::SolverUtils::EquationSystem::m_session.

Referenced by v_DoSolve().

◆ v_DoInitialise()

void Nektar::LaplacePhi::v_DoInitialise ( bool  dumpInitialConditions = true)
overrideprotectedvirtual

Virtual function for initialisation implementation.

By default, nothing needs initialising at the EquationSystem level.

Reimplemented from Nektar::SolverUtils::EquationSystem.

Definition at line 159 of file LaplacePhi.cpp.

160{
161 // Set initial conditions from session file
162 SetInitialConditions(0.0, false, 0);
163}
SOLVER_UTILS_EXPORT void SetInitialConditions(NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
Initialise the data in the dependent fields.

References Nektar::SolverUtils::EquationSystem::SetInitialConditions().

◆ v_DoSolve()

void Nektar::LaplacePhi::v_DoSolve ( void  )
overrideprotectedvirtual

Virtual function for solve implementation.

Reimplemented from Nektar::SolverUtils::EquationSystem.

Definition at line 165 of file LaplacePhi.cpp.

166{
167 for (int i = 0; i < m_fields.size(); ++i)
168 {
170 Array<OneD, NekDouble> forcing(m_fields[i]->GetTotPoints(), 0.);
171 m_fields[i]->HelmSolve(forcing, m_fields[i]->UpdateCoeffs(), m_factors);
172 m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(),
173 m_fields[i]->UpdatePhys());
174 m_fields[i]->SetPhysState(true);
175 }
177}
void setUserDefinedBC(int n)
SOLVER_UTILS_EXPORT int GetTotPoints()

References calculateAddedMass(), Nektar::SolverUtils::EquationSystem::GetTotPoints(), m_factors, Nektar::SolverUtils::EquationSystem::m_fields, and setUserDefinedBC().

◆ v_GenerateSummary()

void Nektar::LaplacePhi::v_GenerateSummary ( SolverUtils::SummaryList l)
overrideprotectedvirtual

Virtual function for generating summary information.

Reimplemented from Nektar::SolverUtils::EquationSystem.

Definition at line 152 of file LaplacePhi.cpp.

153{
155 SolverUtils::AddSummaryItem(s, "Lambda",
157}
SOLVER_UTILS_EXPORT void SessionSummary(SummaryList &vSummary)
Write out a session summary.
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

References Nektar::SolverUtils::AddSummaryItem(), Nektar::StdRegions::eFactorLambda, m_factors, and Nektar::SolverUtils::EquationSystem::SessionSummary().

◆ v_GetSystemSingularChecks()

Array< OneD, bool > Nektar::LaplacePhi::v_GetSystemSingularChecks ( )
overrideprivatevirtual

Reimplemented from Nektar::SolverUtils::EquationSystem.

Definition at line 179 of file LaplacePhi.cpp.

180{
181 return Array<OneD, bool>(m_session->GetVariables().size(), true);
182}

References Nektar::SolverUtils::EquationSystem::m_session.

◆ v_InitObject()

void Nektar::LaplacePhi::v_InitObject ( bool  DeclareFeld = true)
overrideprotectedvirtual

Initialisation object for EquationSystem.

Continuous field

Setting up the normals

Setting up the normals

Reimplemented from Nektar::SolverUtils::EquationSystem.

Definition at line 53 of file LaplacePhi.cpp.

54{
55 EquationSystem::v_InitObject(DeclareFields);
56 int ndim = m_fields[0]->GetCoordim(0);
57 // number of potential fields
58 int bndsize = (ndim == 2) ? 4 : 7;
59 // check variable name phi[0-bndsize]
60 for (const auto &it : m_session->GetVariables())
61 {
62 if (it.size() != 4 || it[0] != 'p' || it[1] != 'h' || it[2] != 'i' ||
63 it[3] < '0' || it[4] >= bndsize + '0')
64 {
66 "Error: incorrect variable name '" + it +
67 "', should be phi[0-3] or phi[0-6].");
68 }
69 }
70 // number of physics points for each field
71 int numpts = 0;
72 Array<OneD, const SpatialDomains::BoundaryConditionShPtr> m_BndConds =
73 m_fields[0]->GetBndConditions();
74 Array<OneD, MultiRegions::ExpListSharedPtr> m_BndExp =
75 m_fields[0]->GetBndCondExpansions();
76 for (int i = 0; i < m_BndConds.size(); ++i)
77 {
78 if (boost::iequals(m_BndConds[i]->GetUserDefined(), "MOVE"))
79 {
80 numpts += m_BndExp[i]->GetTotPoints();
81 }
82 }
83 // load pivot point
84 m_pivot = Array<OneD, NekDouble>(3, 0.);
85 if (m_session->DefinesParameter("Pivotx"))
86 {
87 m_pivot[0] = m_session->GetParameter("Pivotx");
88 }
89 if (m_session->DefinesParameter("Pivoty"))
90 {
91 m_pivot[1] = m_session->GetParameter("Pivoty");
92 }
93 if (m_session->DefinesParameter("Pivotz"))
94 {
95 m_pivot[2] = m_session->GetParameter("Pivotz");
96 }
97 // allocate storage
98 m_boundValues = Array<OneD, Array<OneD, NekDouble>>(bndsize);
99 for (int i = 0; i < bndsize; ++i)
100 {
101 m_boundValues[i] = Array<OneD, NekDouble>(numpts, 0.);
102 }
103 // calculate values
104 int offset = 0;
105 for (int i = 0; i < m_BndConds.size(); ++i)
106 {
107 if (boost::iequals(m_BndConds[i]->GetUserDefined(), "MOVE"))
108 {
109 int npts = m_BndExp[i]->GetTotPoints();
110 Array<OneD, Array<OneD, NekDouble>> n(3);
111 Array<OneD, Array<OneD, NekDouble>> x(3);
112 for (int j = 0; j < 3; ++j)
113 {
114 n[j] = Array<OneD, NekDouble>(npts, 0.);
115 x[j] = Array<OneD, NekDouble>(npts, 0.);
116 }
117 m_BndExp[i]->GetNormals(n);
118 m_BndExp[i]->GetCoords(x[0], x[1], x[2]);
119 for (int j = 0; j < 3; ++j)
120 {
121 Vmath::Sadd(npts, -m_pivot[j], x[j], 1, x[j], 1);
122 }
123 Array<OneD, NekDouble> atmp;
124 for (int j = 0; j < ndim; ++j)
125 {
126 Vmath::Smul(npts, -1., n[j], 1,
127 atmp = m_boundValues[j] + offset, 1);
128 }
129 atmp = m_boundValues[bndsize - 2] + offset;
130 Vmath::Vvtvvtm(npts, &x[1][0], 1, &n[0][0], 1, &x[0][0], 1,
131 &n[1][0], 1, &atmp[0], 1);
132 if (ndim == 3)
133 {
134 atmp = m_boundValues[3] + offset;
135 Vmath::Vvtvvtm(npts, &x[2][0], 1, &n[1][0], 1, &x[1][0], 1,
136 &n[2][0], 1, &atmp[0], 1);
137 atmp = m_boundValues[4] + offset;
138 Vmath::Vvtvvtm(npts, &x[0][0], 1, &n[2][0], 1, &x[2][0], 1,
139 &n[0][0], 1, &atmp[0], 1);
140 }
141 atmp = m_boundValues[bndsize - 1] + offset;
142 for (int j = 0; j < 2; ++j)
143 {
144 Vmath::Vvtvp(npts, x[j], 1, n[j], 1, atmp, 1, atmp, 1);
145 }
146 Vmath::Neg(npts, atmp, 1);
147 offset += npts;
148 }
149 }
150}
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Array< OneD, NekDouble > m_pivot
Definition LaplacePhi.h:82
virtual SOLVER_UTILS_EXPORT void v_InitObject(bool DeclareFeld=true)
Initialisation object for EquationSystem.
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition Vmath.hpp:292
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition Vmath.hpp:366
void Vvtvvtm(int n, const T *v, int incv, const T *w, int incw, const T *x, int incx, const T *y, int incy, T *z, int incz)
vvtvvtm (vector times vector minus vector times vector):
Definition Vmath.hpp:456
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 Sadd(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Add vector y = alpha + x.
Definition Vmath.hpp:194

References Nektar::ErrorUtil::efatal, m_boundValues, Nektar::SolverUtils::EquationSystem::m_fields, m_pivot, Nektar::SolverUtils::EquationSystem::m_session, Vmath::Neg(), NEKERROR, Vmath::Sadd(), Vmath::Smul(), Nektar::SolverUtils::EquationSystem::v_InitObject(), Vmath::Vvtvp(), and Vmath::Vvtvvtm().

Friends And Related Symbol Documentation

◆ MemoryManager< LaplacePhi >

friend class MemoryManager< LaplacePhi >
friend

Class may only be instantiated through the MemoryManager.

Definition at line 1 of file LaplacePhi.h.

Member Data Documentation

◆ className

std::string Nektar::LaplacePhi::className
static
Initial value:
=
static EquationSystemSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Creates an instance of this class.
Definition LaplacePhi.h:51
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
EquationSystemFactory & GetEquationSystemFactory()

Name of class.

Definition at line 62 of file LaplacePhi.h.

◆ m_boundValues

Array<OneD, Array<OneD, NekDouble> > Nektar::LaplacePhi::m_boundValues
private

Definition at line 81 of file LaplacePhi.h.

Referenced by calculateAddedMass(), setUserDefinedBC(), and v_InitObject().

◆ m_factors

StdRegions::ConstFactorMap Nektar::LaplacePhi::m_factors
protected

Definition at line 65 of file LaplacePhi.h.

Referenced by LaplacePhi(), v_DoSolve(), and v_GenerateSummary().

◆ m_pivot

Array<OneD, NekDouble> Nektar::LaplacePhi::m_pivot
private

Definition at line 82 of file LaplacePhi.h.

Referenced by v_InitObject().