Nektar++
Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes | Friends | List of all members
Nektar::SolverUtils::StatLagrangianPoints Class Reference

#include <FilterLagrangianPoints.h>

Inheritance diagram for Nektar::SolverUtils::StatLagrangianPoints:
[legend]

Public Member Functions

 StatLagrangianPoints (int rank, int dim, int intOrder, int idOffset, NekDouble dt, const std::vector< int > &Np, const std::vector< NekDouble > &Box, std::vector< std::string > extraVars)
 
 ~StatLagrangianPoints () override
 
- Public Member Functions inherited from Nektar::SolverUtils::StationaryPoints
virtual ~StationaryPoints ()=default
 
SOLVER_UTILS_EXPORT void OutputData (std::string filename, bool verbose, std::map< std::string, NekDouble > &params)
 
void ReSize (int Np)
 
void AssignPoint (int id, int pid, const Array< OneD, NekDouble > &gcoords)
 
void GetCoordsByPID (int pid, Array< OneD, NekDouble > &gcoords)
 
void SetCoordsByPID (int pid, const Array< OneD, NekDouble > &gcoords)
 
void GetPhysicsByPID (int pid, Array< OneD, NekDouble > &data)
 
void SetPhysicsByPID (int pid, const Array< OneD, NekDouble > &data)
 
void GetCoords (int pid, Array< OneD, NekDouble > &gcoords)
 
void SetCoords (int pid, const Array< OneD, NekDouble > &gcoords)
 
void GetPhysics (int pid, Array< OneD, NekDouble > &data)
 
void SetPhysics (int pid, const Array< OneD, NekDouble > &data)
 
void TimeAdvance (int order)
 
int GetTotPoints ()
 
int GetDim ()
 
int LocalToGlobal (int id)
 
int GlobalToLocal (int id)
 

Static Public Member Functions

static StationaryPointsSharedPtr create (int rank, int dim, int intOrder, int idOffset, NekDouble dt, const std::vector< int > &Np, const std::vector< NekDouble > &Box, std::vector< std::string > extraVars)
 Creates an instance of this class. More...
 

Protected Member Functions

void v_OutputData (std::string filename, bool verbose, std::map< std::string, NekDouble > &params) override
 
void v_TimeAdvance (int order) override
 
void v_GetCoords (int pid, Array< OneD, NekDouble > &gcoords) override
 
void v_SetCoords (int pid, const Array< OneD, NekDouble > &gcoords) override
 
void v_GetPhysics (int pid, Array< OneD, NekDouble > &data) override
 
void v_SetPhysics (int pid, const Array< OneD, NekDouble > &data) override
 
void v_ReSize (int Np) override
 
void v_AssignPoint (int id, int pid, const Array< OneD, NekDouble > &gcoords) override
 
virtual void v_OutputData (std::string filename, bool verbose, std::map< std::string, NekDouble > &params)=0
 
virtual void v_TimeAdvance (int order)=0
 
virtual void v_GetCoords (int pid, Array< OneD, NekDouble > &gcoords)=0
 
virtual void v_SetCoords (int pid, const Array< OneD, NekDouble > &gcoords)=0
 
virtual void v_GetPhysics (int pid, Array< OneD, NekDouble > &data)=0
 
virtual void v_SetPhysics (int pid, const Array< OneD, NekDouble > &data)=0
 
virtual void v_ReSize (int Np)=0
 
virtual void v_AssignPoint (int id, int pid, const Array< OneD, NekDouble > &gcoords)
 

Protected Attributes

Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_coords
 
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_velocity
 
Array< OneD, Array< OneD, NekDouble > > m_extraPhysics
 
std::vector< std::string > m_extraPhysVars
 
NekDouble m_dt
 
int m_idOffset
 
std::vector< int > m_N
 
int m_intOrder
 
- Protected Attributes inherited from Nektar::SolverUtils::StationaryPoints
int m_dim
 
int m_totPts
 
std::map< int, int > m_localIDToGlobal
 
std::map< int, int > m_globalIDToLocal
 

Friends

class MemoryManager< StatLagrangianPoints >
 

Detailed Description

Definition at line 49 of file FilterLagrangianPoints.h.

Constructor & Destructor Documentation

◆ StatLagrangianPoints()

Nektar::SolverUtils::StatLagrangianPoints::StatLagrangianPoints ( int  rank,
int  dim,
int  intOrder,
int  idOffset,
NekDouble  dt,
const std::vector< int > &  Np,
const std::vector< NekDouble > &  Box,
std::vector< std::string >  extraVars 
)

Definition at line 366 of file FilterLagrangianPoints.cpp.

372{
373 m_idOffset = idOffset;
374 m_dim = dim;
375 if (intOrder < 0)
376 {
377 m_intOrder = -intOrder;
378 m_dt = -dt;
379 }
380 else
381 {
382 m_intOrder = intOrder;
383 m_dt = dt;
384 }
385 if (Np.size() < m_dim || Box.size() < m_dim * 2)
386 {
387 throw ErrorUtil::NekError("Not enough input for StationaryPoints.");
388 }
389 Array<OneD, NekDouble> delta(3, 1.);
390 m_N.resize(3, 1);
391 m_totPts = 1;
392 for (int d = 0; d < m_dim; ++d)
393 {
394 m_N[d] = Np[d];
395 m_totPts *= Np[d];
396 delta[d] = Np[d] < 2 ? 1. : (Box[d * 2 + 1] - Box[d * 2]) / (Np[d] - 1);
397 }
398 m_coords = Array<OneD, Array<OneD, Array<OneD, NekDouble>>>(m_intOrder);
399 m_velocity = Array<OneD, Array<OneD, Array<OneD, NekDouble>>>(m_intOrder);
400 for (int i = 0; i < m_intOrder; ++i)
401 {
402 m_coords[i] = Array<OneD, Array<OneD, NekDouble>>(m_dim);
403 m_velocity[i] = Array<OneD, Array<OneD, NekDouble>>(m_dim);
404 for (int d = 0; d < m_dim; ++d)
405 {
406 m_coords[i][d] = Array<OneD, NekDouble>(m_totPts, 0.);
407 m_velocity[i][d] = Array<OneD, NekDouble>(m_totPts, 0.);
408 }
409 }
410 m_extraPhysics = Array<OneD, Array<OneD, NekDouble>>(extraVars.size());
411 for (size_t i = 0; i < extraVars.size(); ++i)
412 {
413 m_extraPhysVars.push_back(extraVars[i]);
414 m_extraPhysics[i] = Array<OneD, NekDouble>(m_totPts, 0.);
415 }
416 // initialise points
417 int count = 0;
418 for (int k = 0; k < m_N[2]; ++k)
419 {
420 for (int j = 0; j < m_N[1]; ++j)
421 {
422 for (int i = 0; i < m_N[0]; ++i)
423 {
424 if (m_dim > 0)
425 {
426 m_coords[0][0][count] = Box[0] + delta[0] * i;
427 }
428 if (m_dim > 1)
429 {
430 m_coords[0][1][count] = Box[2] + delta[1] * j;
431 }
432 if (m_dim > 2)
433 {
434 m_coords[0][2][count] = Box[4] + delta[2] * k;
435 }
436 m_localIDToGlobal[count] = m_idOffset + count;
437 m_globalIDToLocal[m_idOffset + count] = count;
438 ++count;
439 }
440 }
441 }
442}
Nektar::ErrorUtil::NekError NekError
Array< OneD, Array< OneD, NekDouble > > m_extraPhysics
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_velocity
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_coords
std::vector< double > d(NPUPPER *NPUPPER)

References Nektar::UnitTests::d(), m_coords, Nektar::SolverUtils::StationaryPoints::m_dim, m_dt, m_extraPhysics, m_extraPhysVars, Nektar::SolverUtils::StationaryPoints::m_globalIDToLocal, m_idOffset, m_intOrder, Nektar::SolverUtils::StationaryPoints::m_localIDToGlobal, m_N, Nektar::SolverUtils::StationaryPoints::m_totPts, and m_velocity.

◆ ~StatLagrangianPoints()

Nektar::SolverUtils::StatLagrangianPoints::~StatLagrangianPoints ( )
inlineoverride

Definition at line 72 of file FilterLagrangianPoints.h.

73 {
74 }

Member Function Documentation

◆ create()

static StationaryPointsSharedPtr Nektar::SolverUtils::StatLagrangianPoints::create ( int  rank,
int  dim,
int  intOrder,
int  idOffset,
NekDouble  dt,
const std::vector< int > &  Np,
const std::vector< NekDouble > &  Box,
std::vector< std::string >  extraVars 
)
inlinestatic

Creates an instance of this class.

Definition at line 55 of file FilterLagrangianPoints.h.

60 {
63 rank, dim, intOrder, idOffset, dt, Np, Box, extraVars);
64 return p;
65 }
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
std::shared_ptr< StationaryPoints > StationaryPointsSharedPtr

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), and CellMLToNektar.cellml_metadata::p.

◆ v_AssignPoint()

void Nektar::SolverUtils::StatLagrangianPoints::v_AssignPoint ( int  id,
int  pid,
const Array< OneD, NekDouble > &  gcoords 
)
overrideprotectedvirtual

Reimplemented from Nektar::SolverUtils::StationaryPoints.

Definition at line 356 of file FilterLagrangianPoints.cpp.

358{
359 StationaryPoints::v_AssignPoint(id, pid, gcoords);
360 for (int d = 0; d < m_dim; ++d)
361 {
362 m_coords[0][d][id] = gcoords[d];
363 }
364}
virtual void v_AssignPoint(int id, int pid, const Array< OneD, NekDouble > &gcoords)

References Nektar::UnitTests::d(), m_coords, Nektar::SolverUtils::StationaryPoints::m_dim, and Nektar::SolverUtils::StationaryPoints::v_AssignPoint().

◆ v_GetCoords()

void Nektar::SolverUtils::StatLagrangianPoints::v_GetCoords ( int  pid,
Array< OneD, NekDouble > &  gcoords 
)
overrideprotectedvirtual

Implements Nektar::SolverUtils::StationaryPoints.

Definition at line 51 of file FilterLagrangianPoints.cpp.

52{
53 for (int d = 0; d < m_dim; ++d)
54 {
55 gcoords[d] = m_coords[0][d][id];
56 }
57}

References Nektar::UnitTests::d(), m_coords, and Nektar::SolverUtils::StationaryPoints::m_dim.

◆ v_GetPhysics()

void Nektar::SolverUtils::StatLagrangianPoints::v_GetPhysics ( int  pid,
Array< OneD, NekDouble > &  data 
)
overrideprotectedvirtual

Implements Nektar::SolverUtils::StationaryPoints.

Definition at line 68 of file FilterLagrangianPoints.cpp.

69{
70 for (int d = 0; d < m_dim; ++d)
71 {
72 data[d] = m_velocity[0][d][id];
73 }
74}

References Nektar::UnitTests::d(), Nektar::SolverUtils::StationaryPoints::m_dim, and m_velocity.

◆ v_OutputData()

void Nektar::SolverUtils::StatLagrangianPoints::v_OutputData ( std::string  filename,
bool  verbose,
std::map< std::string, NekDouble > &  params 
)
overrideprotectedvirtual

Implements Nektar::SolverUtils::StationaryPoints.

Definition at line 279 of file FilterLagrangianPoints.cpp.

282{
283 std::vector<std::string> COORDS = {"x", "y", "z"};
284 std::vector<std::string> VELOCI = {"u", "v", "w"};
285 std::vector<std::string> variables;
286 std::vector<Array<OneD, NekDouble>> data;
287 for (int d = 0; d < m_dim; ++d)
288 {
289 variables.push_back(COORDS[d]);
290 }
291 for (int d = 0; d < m_dim; ++d)
292 {
293 variables.push_back(VELOCI[d]);
294 }
295 if (params.size() > 0)
296 {
297 for (int d = 0; d < m_dim; ++d)
298 {
299 Array<OneD, NekDouble> x(m_totPts, 0.);
300 Vmath::Sadd(m_totPts, params[COORDS[d]], m_coords[0][d], 1, x, 1);
301 data.push_back(x);
302 }
303 for (int d = 0; d < m_dim; ++d)
304 {
305 Array<OneD, NekDouble> x(m_totPts, 0.);
306 Vmath::Sadd(m_totPts, params[VELOCI[d]], m_velocity[0][d], 1, x, 1);
307 data.push_back(x);
308 }
309 }
310 else
311 {
312 for (int d = 0; d < m_dim; ++d)
313 {
314 data.push_back(m_coords[0][d]);
315 }
316 for (int d = 0; d < m_dim; ++d)
317 {
318 data.push_back(m_velocity[0][d]);
319 }
320 }
321 for (size_t d = 0; d < m_extraPhysics.size(); ++d)
322 {
323 data.push_back(m_extraPhysics[d]);
324 variables.push_back(m_extraPhysVars[d]);
325 }
326 OutputTec360_binary(filename, variables, m_N, data, 1);
327 if (verbose)
328 {
329 int Np = m_coords[0][0].size();
330 for (int d = 0; d < m_dim; ++d)
331 {
332 NekDouble value = 0.;
333 for (int i = 0; i < Np; ++i)
334 {
335 value += m_coords[0][d][i] * m_coords[0][d][i];
336 }
337 value = sqrt(value / Np);
338 cout << "L 2 error (variable L" << COORDS[d] << ") : " << value
339 << endl;
340 }
341 for (int d = 0; d < m_dim; ++d)
342 {
343 NekDouble value = 0.;
344 for (int i = 0; i < Np; ++i)
345 {
346 value += m_velocity[0][d][i] * m_velocity[0][d][i];
347 }
348 value = sqrt(value / Np);
349 cout << "L 2 error (variable L" << VELOCI[d] << ") : " << value
350 << endl;
351 }
352 }
353 cout << "Write file " << filename << endl;
354}
static int OutputTec360_binary(const std::string filename, const std::vector< std::string > &variables, const std::vector< int > &rawN, std::vector< Array< OneD, NekDouble > > &data, int isdouble)
double NekDouble
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
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294

References Nektar::UnitTests::d(), m_coords, Nektar::SolverUtils::StationaryPoints::m_dim, m_extraPhysics, m_extraPhysVars, m_N, Nektar::SolverUtils::StationaryPoints::m_totPts, m_velocity, Nektar::SolverUtils::OutputTec360_binary(), Vmath::Sadd(), and tinysimd::sqrt().

◆ v_ReSize()

void Nektar::SolverUtils::StatLagrangianPoints::v_ReSize ( int  Np)
overrideprotectedvirtual

Implements Nektar::SolverUtils::StationaryPoints.

Definition at line 135 of file FilterLagrangianPoints.cpp.

136{
137 if (Np > m_totPts)
138 {
139 for (int i = 0; i < m_coords.size(); ++i)
140 {
141 for (int d = 0; d < m_dim; ++d)
142 {
143 m_coords[i][d] = Array<OneD, NekDouble>(Np);
144 m_velocity[i][d] = Array<OneD, NekDouble>(Np);
145 }
146 }
147 }
148 m_totPts = Np;
149 m_localIDToGlobal.clear();
150 m_globalIDToLocal.clear();
151}

References Nektar::UnitTests::d(), m_coords, Nektar::SolverUtils::StationaryPoints::m_dim, Nektar::SolverUtils::StationaryPoints::m_globalIDToLocal, Nektar::SolverUtils::StationaryPoints::m_localIDToGlobal, Nektar::SolverUtils::StationaryPoints::m_totPts, and m_velocity.

◆ v_SetCoords()

void Nektar::SolverUtils::StatLagrangianPoints::v_SetCoords ( int  pid,
const Array< OneD, NekDouble > &  gcoords 
)
overrideprotectedvirtual

Implements Nektar::SolverUtils::StationaryPoints.

Definition at line 59 of file FilterLagrangianPoints.cpp.

61{
62 for (int d = 0; d < m_dim; ++d)
63 {
64 m_coords[0][d][id] = gcoords[d];
65 }
66}

References Nektar::UnitTests::d(), m_coords, and Nektar::SolverUtils::StationaryPoints::m_dim.

◆ v_SetPhysics()

void Nektar::SolverUtils::StatLagrangianPoints::v_SetPhysics ( int  pid,
const Array< OneD, NekDouble > &  data 
)
overrideprotectedvirtual

Implements Nektar::SolverUtils::StationaryPoints.

Definition at line 76 of file FilterLagrangianPoints.cpp.

78{
79 for (int d = 0; d < m_dim; ++d)
80 {
81 m_velocity[0][d][id] = data[d];
82 }
83 for (size_t d = 0; d < m_extraPhysVars.size(); ++d)
84 {
85 m_extraPhysics[d][id] = data[d + m_dim];
86 }
87}

References Nektar::UnitTests::d(), Nektar::SolverUtils::StationaryPoints::m_dim, m_extraPhysics, m_extraPhysVars, and m_velocity.

◆ v_TimeAdvance()

void Nektar::SolverUtils::StatLagrangianPoints::v_TimeAdvance ( int  order)
overrideprotectedvirtual

Implements Nektar::SolverUtils::StationaryPoints.

Definition at line 104 of file FilterLagrangianPoints.cpp.

105{
106 order = min(order, m_intOrder);
107 if (order == 1)
108 {
109 for (int d = 0; d < m_dim; ++d)
110 {
112 m_coords[0][d], 1);
113 }
115 }
116 else if (order == 2)
117 {
118 for (int d = 0; d < m_dim; ++d)
119 {
120 Array<OneD, NekDouble> res = m_coords[1][d];
121 Vmath::Svtvm(m_totPts, 3., m_velocity[0][d], 1, m_velocity[1][d], 1,
122 res, 1);
123 Vmath::Svtvp(m_totPts, 0.5 * m_dt, res, 1, m_coords[0][d], 1, res,
124 1);
125 }
128 }
129 else
130 {
131 ASSERTL0(false, "Integration order not supported.");
132 }
133}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
static void RollOver(Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &data)
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 Svtvm(int n, const T alpha, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Svtvm (scalar times vector minus vector): z = alpha*x - y.
Definition: Vmath.hpp:424

References ASSERTL0, Nektar::UnitTests::d(), m_coords, Nektar::SolverUtils::StationaryPoints::m_dim, m_dt, m_intOrder, Nektar::SolverUtils::StationaryPoints::m_totPts, m_velocity, Nektar::SolverUtils::RollOver(), Vmath::Svtvm(), and Vmath::Svtvp().

Friends And Related Function Documentation

◆ MemoryManager< StatLagrangianPoints >

friend class MemoryManager< StatLagrangianPoints >
friend

Definition at line 1 of file FilterLagrangianPoints.h.

Member Data Documentation

◆ m_coords

Array<OneD, Array<OneD, Array<OneD, NekDouble> > > Nektar::SolverUtils::StatLagrangianPoints::m_coords
protected

◆ m_dt

NekDouble Nektar::SolverUtils::StatLagrangianPoints::m_dt
protected

Definition at line 91 of file FilterLagrangianPoints.h.

Referenced by StatLagrangianPoints(), and v_TimeAdvance().

◆ m_extraPhysics

Array<OneD, Array<OneD, NekDouble> > Nektar::SolverUtils::StatLagrangianPoints::m_extraPhysics
protected

Definition at line 89 of file FilterLagrangianPoints.h.

Referenced by StatLagrangianPoints(), v_OutputData(), and v_SetPhysics().

◆ m_extraPhysVars

std::vector<std::string> Nektar::SolverUtils::StatLagrangianPoints::m_extraPhysVars
protected

Definition at line 90 of file FilterLagrangianPoints.h.

Referenced by StatLagrangianPoints(), v_OutputData(), and v_SetPhysics().

◆ m_idOffset

int Nektar::SolverUtils::StatLagrangianPoints::m_idOffset
protected

Definition at line 92 of file FilterLagrangianPoints.h.

Referenced by StatLagrangianPoints().

◆ m_intOrder

int Nektar::SolverUtils::StatLagrangianPoints::m_intOrder
protected

Definition at line 94 of file FilterLagrangianPoints.h.

Referenced by StatLagrangianPoints(), and v_TimeAdvance().

◆ m_N

std::vector<int> Nektar::SolverUtils::StatLagrangianPoints::m_N
protected

Definition at line 93 of file FilterLagrangianPoints.h.

Referenced by StatLagrangianPoints(), and v_OutputData().

◆ m_velocity

Array<OneD, Array<OneD, Array<OneD, NekDouble> > > Nektar::SolverUtils::StatLagrangianPoints::m_velocity
protected