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

#include <EvaluatePoints.h>

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

Public Member Functions

SOLVER_UTILS_EXPORT EvaluatePoints ()
 
virtual ~EvaluatePoints ()
 
SOLVER_UTILS_EXPORT void PartitionMobilePts (const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields)
 
SOLVER_UTILS_EXPORT void EvaluateMobilePhys (const MultiRegions::ExpListSharedPtr &pField, std::vector< Array< OneD, NekDouble > > &PhysicsData, NekDouble time)
 
SOLVER_UTILS_EXPORT void PassMobilePhysToStatic (std::map< int, std::set< int > > &callbackUpdateMobCoords)
 
SOLVER_UTILS_EXPORT void PassStaticCoordsToMobile (std::map< int, std::set< int > > &callbackUpdateMobCoords)
 
SOLVER_UTILS_EXPORT void CopyStaticPtsToMobile ()
 
SOLVER_UTILS_EXPORT void CopyMobilePtsToStatic ()
 
SOLVER_UTILS_EXPORT void Initialise (const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time, std::vector< std::string > &defaultValues)
 
SOLVER_UTILS_EXPORT void SetUpCommInfo ()
 

Protected Member Functions

void UpdateMobileCoords (std::map< int, Array< OneD, NekDouble > > &globalCoords)
 
void GatherMobilePhysics (std::map< int, Array< OneD, NekDouble > > &revData)
 
void PartitionLocalPoints (const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, std::set< int > &notfound)
 
void PartitionExchangeNonlocalPoints (const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, std::set< int > &notfound)
 
void SyncColumnComm ()
 
NekDouble GetDefaultValue (const int i, const Array< OneD, const NekDouble > x, const NekDouble time)
 
void Pack2Int (const int &a, const int &b, double &d)
 
void unPack2Int (int &a, int &b, const double &d)
 
void Pack3Int (const int &a, const int &b, const int &c, double &d)
 
void unPack3Int (int &a, int &b, int &c, const double &d)
 
void Pack2Short (const int &a, const int &b, int &c)
 
void unPack2Short (int &a, int &b, const int &c)
 
virtual void v_ModifyVelocity (Array< OneD, NekDouble > gcoords, NekDouble time, Array< OneD, NekDouble > vel)
 

Protected Attributes

int m_spacedim
 
int m_rank
 
int m_rowRank
 
int m_columnRank
 
bool m_isHomogeneous1D
 
int m_nPhysics
 
LibUtilities::CommSharedPtr m_comm
 
StationaryPointsSharedPtr m_staticPts
 
std::map< int, MobilePointSharedPtrm_mobilePts
 
std::map< int, LibUtilities::EquationSharedPtrm_defaultValueFunction
 
std::map< int, int > m_recvMobInfo
 
std::map< int, int > m_recvStatInfo
 

Detailed Description

Definition at line 210 of file EvaluatePoints.h.

Constructor & Destructor Documentation

◆ EvaluatePoints()

Nektar::SolverUtils::EvaluatePoints::EvaluatePoints ( )

Definition at line 104 of file EvaluatePoints.cpp.

105{
106}

◆ ~EvaluatePoints()

Nektar::SolverUtils::EvaluatePoints::~EvaluatePoints ( )
virtual

Definition at line 111 of file EvaluatePoints.cpp.

112{
113}

Member Function Documentation

◆ CopyMobilePtsToStatic()

void Nektar::SolverUtils::EvaluatePoints::CopyMobilePtsToStatic ( )

Definition at line 737 of file EvaluatePoints.cpp.

738{
740 "CopyStaticPtsToMobile shoule only be called by ColumnRank zero");
741 m_staticPts->ReSize(m_mobilePts.size());
742 int id = 0;
743 for (const auto &p : m_mobilePts)
744 {
745 m_staticPts->AssignPoint(id, p.first, p.second->GetGlobalCoords());
746 }
747}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
StationaryPointsSharedPtr m_staticPts
std::map< int, MobilePointSharedPtr > m_mobilePts

References ASSERTL0, m_columnRank, m_isHomogeneous1D, m_mobilePts, m_staticPts, and CellMLToNektar.cellml_metadata::p.

◆ CopyStaticPtsToMobile()

void Nektar::SolverUtils::EvaluatePoints::CopyStaticPtsToMobile ( )

Definition at line 722 of file EvaluatePoints.cpp.

723{
725 "CopyStaticPtsToMobile shoule only be called by ColumnRank zero");
726 m_mobilePts.clear();
727 for (int p = 0; p < m_staticPts->GetTotPoints(); ++p)
728 {
729 Array<OneD, NekDouble> data(m_spacedim, 0.);
730 m_staticPts->GetCoords(p, data);
731 m_mobilePts[m_staticPts->LocalToGlobal(p)] =
733 data);
734 }
735}
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), ASSERTL0, m_columnRank, m_isHomogeneous1D, m_mobilePts, m_rank, m_spacedim, m_staticPts, and CellMLToNektar.cellml_metadata::p.

Referenced by Nektar::SolverUtils::FilterLagrangianPoints::v_Initialise().

◆ EvaluateMobilePhys()

void Nektar::SolverUtils::EvaluatePoints::EvaluateMobilePhys ( const MultiRegions::ExpListSharedPtr pField,
std::vector< Array< OneD, NekDouble > > &  PhysicsData,
NekDouble  time 
)

Definition at line 154 of file EvaluatePoints.cpp.

157{
158 m_nPhysics = PhysicsData.size();
159 int coordim = pField->GetGraph()->GetSpaceDimension();
160 Array<OneD, NekDouble> data(m_mobilePts.size() * m_nPhysics, 0.0);
161 Array<OneD, NekDouble> physvals;
162 Array<OneD, NekDouble> locCoord;
164 {
165 Array<OneD, const unsigned int> planes = pField->GetZIDs();
166 int nPlanes = pField->GetHomogeneousBasis()->GetZ().size();
167 NekDouble lHom = pField->GetHomoLen();
168 int nelmts = pField->GetPlane(0)->GetExpSize();
169 int indx = -1;
170 for (auto &x : m_mobilePts)
171 {
172 ++indx;
173 locCoord = x.second->GetLocalCoords();
174 int expId = x.second->m_eId;
175 if (expId < 0)
176 {
177 continue;
178 }
179 NekDouble BetaT = 2. * M_PI * fmod(locCoord[coordim], lHom) / lHom;
180 std::vector<NekDouble> basis(planes.size(), 0.);
181 for (size_t n = 0; n < planes.size(); ++n)
182 {
183 if (planes[n] == 0)
184 {
185 basis[n] = 1.;
186 }
187 else if (planes[n] == 1)
188 {
189 basis[n] = cos(0.5 * nPlanes * BetaT);
190 }
191 else if (planes[n] % 2 == 0)
192 {
193 NekDouble phase = (planes[n] >> 1) * BetaT;
194 basis[n] = cos(phase);
195 }
196 else
197 {
198 NekDouble phase = (planes[n] >> 1) * BetaT;
199 basis[n] = -sin(phase);
200 }
201 }
203 pField->GetPlane(0)->GetExp(expId);
204 for (int j = 0; j < m_nPhysics; ++j)
205 {
206 NekDouble value = 0.0;
207 for (size_t n = 0; n < planes.size(); ++n)
208 {
209 physvals = PhysicsData[j] +
210 pField->GetPhys_Offset(expId + n * nelmts);
211 NekDouble coeff = elmt->StdPhysEvaluate(locCoord, physvals);
212 value += basis[n] * coeff;
213 }
214 data[indx * m_nPhysics + j] = value;
215 }
216 }
217 m_comm->GetColumnComm()->AllReduce(data, LibUtilities::ReduceSum);
218 }
219 else
220 {
221 for (int j = 0; j < m_nPhysics; ++j)
222 {
223 int indx = -1;
224 for (auto &x : m_mobilePts)
225 {
226 ++indx;
227 locCoord = x.second->GetLocalCoords();
228 int expId = x.second->m_eId;
229 if (expId < 0)
230 {
231 continue;
232 }
233 physvals = PhysicsData[j] + pField->GetPhys_Offset(expId);
234 data[indx * m_nPhysics + j] =
235 pField->GetExp(expId)->StdPhysEvaluate(locCoord, physvals);
236 }
237 }
238 }
239 // copy data
240 int count = 0;
241 for (auto &x : m_mobilePts)
242 {
243 Array<OneD, NekDouble> datum;
244 if (x.second->m_eId >= 0)
245 {
246 datum = data + count * m_nPhysics;
247 }
248 else
249 {
250 datum = Array<OneD, NekDouble>(m_nPhysics, 0.);
251 if (m_columnRank == 0)
252 {
253 for (int j = 0; j < m_nPhysics; ++j)
254 {
255 datum[j] =
256 GetDefaultValue(j, x.second->GetGlobalCoords(), time);
257 }
258 }
259 }
260 v_ModifyVelocity(x.second->GetGlobalCoords(), time, datum);
261 x.second->SetData(m_nPhysics, datum);
262 ++count;
263 }
264}
NekDouble GetDefaultValue(const int i, const Array< OneD, const NekDouble > x, const NekDouble time)
virtual void v_ModifyVelocity(Array< OneD, NekDouble > gcoords, NekDouble time, Array< OneD, NekDouble > vel)
LibUtilities::CommSharedPtr m_comm
std::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:66
double NekDouble

References GetDefaultValue(), m_columnRank, m_comm, m_isHomogeneous1D, m_mobilePts, m_nPhysics, Nektar::LibUtilities::ReduceSum, and v_ModifyVelocity().

Referenced by Nektar::SolverUtils::FilterLagrangianPoints::v_Update().

◆ GatherMobilePhysics()

void Nektar::SolverUtils::EvaluatePoints::GatherMobilePhysics ( std::map< int, Array< OneD, NekDouble > > &  revData)
protected

Definition at line 274 of file EvaluatePoints.cpp.

276{
277 revData.clear();
278 // send/receive data
279 int nPackageSize = m_nPhysics + 1;
280 if (m_columnRank == 0)
281 {
282 std::map<int, std::set<int>> threadToSend;
283 for (const auto &p : m_mobilePts)
284 {
285 threadToSend[p.second->m_sRank].insert(p.first);
286 }
287 for (const auto &thread : threadToSend)
288 {
289 Array<OneD, NekDouble> dataToSend(
290 thread.second.size() * nPackageSize, 0.);
291 Array<OneD, NekDouble> tmp;
292 int offset = 0;
293 for (const int p : thread.second)
294 {
295 Pack3Int(p, m_rank, thread.first, dataToSend[offset]);
296 Vmath::Vcopy(m_nPhysics, m_mobilePts[p]->GetData(), 1,
297 tmp = dataToSend + offset + 1, 1);
298 offset += nPackageSize;
299 }
300 if (thread.first != m_rank)
301 {
302 m_comm->Send(thread.first, dataToSend);
303 }
304 else
305 {
306 revData[thread.first] = dataToSend;
307 }
308 }
309 }
310 for (const auto &p : m_recvMobInfo)
311 {
312 if (p.first != m_rank)
313 {
314 Array<OneD, NekDouble> dataToRecv(p.second * nPackageSize, 0.);
315 m_comm->Recv(p.first, dataToRecv);
316 revData[p.first] = dataToRecv;
317 }
318 }
319}
void Pack3Int(const int &a, const int &b, const int &c, double &d)
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825

References m_columnRank, m_comm, m_mobilePts, m_nPhysics, m_rank, m_recvMobInfo, CellMLToNektar.cellml_metadata::p, Pack3Int(), and Vmath::Vcopy().

Referenced by PassMobilePhysToStatic().

◆ GetDefaultValue()

NekDouble Nektar::SolverUtils::EvaluatePoints::GetDefaultValue ( const int  i,
const Array< OneD, const NekDouble x,
const NekDouble  time 
)
protected

Definition at line 673 of file EvaluatePoints.cpp.

676{
677 NekDouble z = m_spacedim == 2 ? 0 : x[2];
679 {
680 return 0.;
681 }
682 else
683 {
684 return m_defaultValueFunction[i]->Evaluate(x[0], x[1], z, time);
685 }
686}
std::map< int, LibUtilities::EquationSharedPtr > m_defaultValueFunction
std::vector< double > z(NPUPPER)

References m_defaultValueFunction, m_spacedim, and Nektar::UnitTests::z().

Referenced by EvaluateMobilePhys().

◆ Initialise()

void Nektar::SolverUtils::EvaluatePoints::Initialise ( const Array< OneD, const MultiRegions::ExpListSharedPtr > &  pFields,
const NekDouble time,
std::vector< std::string > &  defaultValues 
)

Definition at line 118 of file EvaluatePoints.cpp.

122{
123 m_comm = pFields[0]->GetComm();
124 m_rank = m_comm->GetRank();
125 m_columnRank = m_comm->GetColumnComm()->GetRank();
126 m_rowRank = m_comm->GetRowComm()->GetRank();
127 LibUtilities::SessionReaderSharedPtr session = pFields[0]->GetSession();
128 session->MatchSolverInfo("Homogeneous", "1D", m_isHomogeneous1D, false);
129 int coordim = pFields[0]->GetGraph()->GetSpaceDimension();
130 int numHomoDir = m_isHomogeneous1D ? 1 : 0;
131 m_spacedim = coordim + numHomoDir;
132 ASSERTL0(m_rank <= 65535,
133 "EvaluatePoints does not support more than 65536 threads.")
134
135 for (size_t i = 0; i < defaultValues.size(); ++i)
136 {
137 try
138 {
141 session->GetInterpreter(), defaultValues[i]);
142 }
143 catch (const std::runtime_error &)
144 {
146 "Error parsering expression" + defaultValues[i]);
147 }
148 }
149}
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:202
std::shared_ptr< SessionReader > SessionReaderSharedPtr

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), ASSERTL0, Nektar::ErrorUtil::efatal, m_columnRank, m_comm, m_defaultValueFunction, m_isHomogeneous1D, m_rank, m_rowRank, m_spacedim, and NEKERROR.

Referenced by Nektar::SolverUtils::FilterLagrangianPoints::v_Initialise().

◆ Pack2Int()

void Nektar::SolverUtils::EvaluatePoints::Pack2Int ( const int &  a,
const int &  b,
double &  d 
)
protected

Definition at line 55 of file EvaluatePoints.cpp.

56{
57 int *p = (int *)&d;
58 p[0] = a;
59 p[1] = b;
60}
std::vector< double > d(NPUPPER *NPUPPER)

References Nektar::UnitTests::d(), and CellMLToNektar.cellml_metadata::p.

Referenced by PartitionExchangeNonlocalPoints(), and SyncColumnComm().

◆ Pack2Short()

void Nektar::SolverUtils::EvaluatePoints::Pack2Short ( const int &  a,
const int &  b,
int &  c 
)
protected

Definition at line 88 of file EvaluatePoints.cpp.

89{
90 short int *p = (short int *)&c;
91 p[0] = (short int)a;
92 p[1] = (short int)b;
93}

References CellMLToNektar.cellml_metadata::p.

Referenced by SetUpCommInfo().

◆ Pack3Int()

void Nektar::SolverUtils::EvaluatePoints::Pack3Int ( const int &  a,
const int &  b,
const int &  c,
double &  d 
)
protected

Definition at line 69 of file EvaluatePoints.cpp.

71{
72 int *p = (int *)&d;
73 unsigned short int *q = (unsigned short int *)&d;
74 p[0] = a;
75 q[2] = (unsigned short int)b;
76 q[3] = (unsigned short int)c;
77}
std::vector< double > q(NPUPPER *NPUPPER)

References Nektar::UnitTests::d(), CellMLToNektar.cellml_metadata::p, and Nektar::UnitTests::q().

Referenced by GatherMobilePhysics(), PartitionExchangeNonlocalPoints(), and PassStaticCoordsToMobile().

◆ PartitionExchangeNonlocalPoints()

void Nektar::SolverUtils::EvaluatePoints::PartitionExchangeNonlocalPoints ( const Array< OneD, const MultiRegions::ExpListSharedPtr > &  pFields,
std::set< int > &  notfound 
)
protected

Definition at line 464 of file EvaluatePoints.cpp.

467{
468 // Determine the unique process responsible for each history point
469 // For points on a partition boundary, must select a single process
470 int nExchange = 0, nTot = 0;
471 Array<OneD, NekDouble> exchangeInfo;
472 std::vector<int> npts(m_comm->GetRowComm()->GetSize(), 0);
473 if (m_columnRank == 0 || !m_isHomogeneous1D)
474 {
475 npts[m_rowRank] = notfound.size();
476 m_comm->GetRowComm()->AllReduce(npts, LibUtilities::ReduceMax);
477 for (auto v : npts)
478 {
479 nTot += v;
480 }
481 }
482 m_comm->AllReduce(nTot, LibUtilities::ReduceMax);
483 if (nTot == 0)
484 {
485 return;
486 }
487 if (m_columnRank == 0 || !m_isHomogeneous1D)
488 {
489 // spread unfound points to all Row threads
490 int nOffset = 0;
491 for (int i = 0; i < m_rowRank; ++i)
492 {
493 nOffset += npts[i];
494 }
495 int nDataPackage = 1 + m_spacedim;
496 int count = nOffset * nDataPackage;
497 Array<OneD, NekDouble> data(nTot * nDataPackage, 0.);
498 for (auto &p : notfound)
499 {
500 Pack2Int(p, m_mobilePts[p]->m_sRank, data[count++]);
501 for (int i = 0; i < m_spacedim; ++i)
502 {
503 data[count++] = m_mobilePts[p]->GetGlobalCoords()[i];
504 }
505 }
506 m_comm->GetRowComm()->AllReduce(data, LibUtilities::ReduceSum);
507 // global search
508 std::vector<int> eIds(nTot, -1);
509 std::vector<int> procIds(nTot, -1);
510 Array<OneD, NekDouble> localCoords(nTot * m_spacedim, 0.);
511 for (int p = 0; p < nTot; ++p)
512 {
513 if (nOffset <= p && p < nOffset + npts[m_rowRank])
514 {
515 continue;
516 }
517 Array<OneD, NekDouble> gloCoords = data + p * nDataPackage + 1;
518 Array<OneD, NekDouble> locCoords = localCoords + p * m_spacedim;
519 // Determine the expansion and local coordinates
521 {
522 eIds[p] =
523 pFields[0]->GetPlane(0)->GetExpIndex(gloCoords, locCoords);
524 }
525 else
526 {
527 eIds[p] = pFields[0]->GetExpIndex(gloCoords, locCoords);
528 }
529 if (eIds[p] >= 0)
530 {
531 procIds[p] = m_rowRank;
532 }
533 }
534 m_comm->GetRowComm()->AllReduce(procIds, LibUtilities::ReduceMax);
535 for (int p = 0; p < nTot; ++p)
536 {
537 if (procIds[p] >= 0)
538 {
539 ++nExchange;
540 if (procIds[p] == m_rowRank)
541 {
542 int pId, sRank;
543 unPack2Int(pId, sRank, data[p * nDataPackage]);
544 m_mobilePts[pId] =
546 m_spacedim, sRank, data + p * nDataPackage + 1);
547 m_mobilePts[pId]->SetLocalCoords(
548 eIds[p], localCoords + p * m_spacedim);
549 m_recvStatInfo[sRank] += 1;
550 }
551 }
552 }
553 for (int p = nOffset; p < nOffset + npts[m_rowRank]; ++p)
554 {
555 if (procIds[p] >= 0)
556 {
557 int pId, sRank;
558 unPack2Int(pId, sRank, data[p * nDataPackage]);
559 m_mobilePts.erase(pId);
560 if (m_recvStatInfo[sRank] == 1)
561 {
562 m_recvStatInfo.erase(sRank);
563 }
564 else
565 {
566 m_recvStatInfo[sRank] -= 1;
567 }
568 }
569 }
570 // build exchangeInfo [pid, from thread rank, to thread rank]
571 int accumulate = 0;
572 exchangeInfo = Array<OneD, NekDouble>(nExchange, 0.);
573 for (int t = 0, p = 0; t < npts.size(); ++t)
574 {
575 for (int i = 0; i < npts[t]; ++i, ++p)
576 {
577 if (procIds[p] >= 0)
578 {
579 int pId, sRank;
580 unPack2Int(pId, sRank, data[p * nDataPackage]);
581 Pack3Int(sRank, t, procIds[p], exchangeInfo[accumulate++]);
582 }
583 }
584 }
585 }
586 m_comm->AllReduce(nExchange, LibUtilities::ReduceMax);
587 if (nExchange == 0)
588 {
589 return;
590 }
591 // sync column communication
592 if (m_comm->GetColumnComm()->GetSize() > 1)
593 {
594 if (m_columnRank > 0)
595 {
596 exchangeInfo = Array<OneD, NekDouble>(nExchange, 0.);
597 }
598 m_comm->GetColumnComm()->AllReduce(exchangeInfo,
600 }
601 // update m_recvMobInfo
602 for (int p = 0; p < nExchange; ++p)
603 {
604 int sRank, fromRank, toRank;
605 unPack3Int(sRank, fromRank, toRank, exchangeInfo[p]);
606 if (sRank == m_rank)
607 {
608 if (m_recvMobInfo[fromRank] == 1)
609 {
610 m_recvMobInfo.erase(fromRank);
611 }
612 else
613 {
614 m_recvMobInfo[fromRank] -= 1;
615 }
616 m_recvMobInfo[toRank] += 1;
617 }
618 }
619}
void Pack2Int(const int &a, const int &b, double &d)
void unPack2Int(int &a, int &b, const double &d)
void unPack3Int(int &a, int &b, int &c, const double &d)

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), m_columnRank, m_comm, m_isHomogeneous1D, m_mobilePts, m_rank, m_recvMobInfo, m_recvStatInfo, m_rowRank, m_spacedim, CellMLToNektar.cellml_metadata::p, Pack2Int(), Pack3Int(), Nektar::LibUtilities::ReduceMax, Nektar::LibUtilities::ReduceSum, unPack2Int(), and unPack3Int().

Referenced by PartitionMobilePts().

◆ PartitionLocalPoints()

void Nektar::SolverUtils::EvaluatePoints::PartitionLocalPoints ( const Array< OneD, const MultiRegions::ExpListSharedPtr > &  pFields,
std::set< int > &  notfound 
)
protected

Definition at line 429 of file EvaluatePoints.cpp.

432{
434 {
435 return;
436 }
437 // not repeat in local points
438 for (auto &p : m_mobilePts)
439 {
440 Array<OneD, NekDouble> locCoords = p.second->GetLocalCoords();
441 Array<OneD, NekDouble> gloCoords = p.second->GetGlobalCoords();
442 // Determine the expansion and local coordinates
444 {
445 p.second->m_eId = pFields[0]->GetPlane(0)->GetExpIndex(
446 gloCoords, locCoords, 0., false, p.second->m_eId);
447 }
448 else
449 {
450 p.second->m_eId = pFields[0]->GetExpIndex(gloCoords, locCoords, 0.,
451 false, p.second->m_eId);
452 }
453 if (p.second->m_eId < 0)
454 {
455 notfound.insert(p.first);
456 }
457 }
458}

References m_columnRank, m_isHomogeneous1D, m_mobilePts, and CellMLToNektar.cellml_metadata::p.

Referenced by PartitionMobilePts().

◆ PartitionMobilePts()

void Nektar::SolverUtils::EvaluatePoints::PartitionMobilePts ( const Array< OneD, const MultiRegions::ExpListSharedPtr > &  pFields)

Definition at line 381 of file EvaluatePoints.cpp.

383{
384 std::set<int> notfound;
385 PartitionLocalPoints(pFields, notfound);
386 PartitionExchangeNonlocalPoints(pFields, notfound);
388}
void PartitionLocalPoints(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, std::set< int > &notfound)
void PartitionExchangeNonlocalPoints(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, std::set< int > &notfound)

References PartitionExchangeNonlocalPoints(), PartitionLocalPoints(), and SyncColumnComm().

Referenced by Nektar::SolverUtils::FilterLagrangianPoints::v_Update().

◆ PassMobilePhysToStatic()

void Nektar::SolverUtils::EvaluatePoints::PassMobilePhysToStatic ( std::map< int, std::set< int > > &  callbackUpdateMobCoords)

Definition at line 390 of file EvaluatePoints.cpp.

392{
393 std::map<int, Array<OneD, NekDouble>> recvData;
394 GatherMobilePhysics(recvData);
395 // unpack data
396 for (auto const &data : recvData)
397 {
398 int offset = 0;
399 while (offset < data.second.size())
400 {
401 int pid, mobt, statt;
402 unPack3Int(pid, mobt, statt, data.second[offset++]);
403 m_staticPts->SetPhysicsByPID(pid, data.second + offset);
404 callbackUpdateMobCoords[data.first].insert(pid);
405 offset += m_nPhysics;
406 }
407 }
408}
void GatherMobilePhysics(std::map< int, Array< OneD, NekDouble > > &revData)

References GatherMobilePhysics(), m_nPhysics, m_staticPts, and unPack3Int().

Referenced by Nektar::SolverUtils::FilterLagrangianPoints::v_Update().

◆ PassStaticCoordsToMobile()

void Nektar::SolverUtils::EvaluatePoints::PassStaticCoordsToMobile ( std::map< int, std::set< int > > &  callbackUpdateMobCoords)

Definition at line 410 of file EvaluatePoints.cpp.

412{
413 std::map<int, Array<OneD, NekDouble>> globalCoords;
414 for (auto const &thread : callbackUpdateMobCoords)
415 {
416 int nDataSize = m_spacedim + 1, offset = 0;
417 Array<OneD, NekDouble> data(thread.second.size() * nDataSize), tmp;
418 for (const auto &p : thread.second)
419 {
420 Pack3Int(p, thread.first, m_rank, data[offset++]);
421 m_staticPts->GetCoordsByPID(p, tmp = data + offset);
422 offset += m_spacedim;
423 }
424 globalCoords[thread.first] = data;
425 }
426 UpdateMobileCoords(globalCoords);
427}
void UpdateMobileCoords(std::map< int, Array< OneD, NekDouble > > &globalCoords)

References m_rank, m_spacedim, m_staticPts, CellMLToNektar.cellml_metadata::p, Pack3Int(), and UpdateMobileCoords().

Referenced by Nektar::SolverUtils::FilterLagrangianPoints::v_Update().

◆ SetUpCommInfo()

void Nektar::SolverUtils::EvaluatePoints::SetUpCommInfo ( )

Definition at line 321 of file EvaluatePoints.cpp.

322{
323 // setup communication pair;
324 int nRows = m_comm->GetRowComm()->GetSize();
325 int nColumn = m_comm->GetColumnComm()->GetSize();
326 //(from short/int, to short/int, size/int)
327 m_recvStatInfo.clear();
328 int totPairs = 0;
329 Array<OneD, int> data;
330 if (m_columnRank == 0)
331 {
332 for (const auto &p : m_mobilePts)
333 {
334 m_recvStatInfo[p.second->m_sRank] += 1;
335 }
336 Array<OneD, int> npairs(nRows, 0);
337 npairs[m_rowRank] = m_recvStatInfo.size();
338 m_comm->GetRowComm()->AllReduce(npairs, LibUtilities::ReduceMax);
339 int offset = 0;
340 for (int t = 0; t < m_rowRank; ++t)
341 {
342 offset += npairs[t];
343 }
344 totPairs = offset;
345 for (int t = m_rowRank; t < nRows; ++t)
346 {
347 totPairs += npairs[t];
348 }
349 data = Array<OneD, int>(totPairs * 2, 0);
350 offset *= 2;
351 for (const auto &p : m_recvStatInfo)
352 {
353 Pack2Short((short int)m_rowRank, (short int)p.first,
354 data[offset++]);
355 data[offset++] = p.second;
356 }
357 m_comm->GetRowComm()->AllReduce(data, LibUtilities::ReduceSum);
358 }
359 // set up to-receive mobile points
360 if (nColumn > 1)
361 {
362 m_comm->GetColumnComm()->AllReduce(totPairs, LibUtilities::ReduceMax);
363 if (m_columnRank > 0)
364 {
365 data = Array<OneD, int>(totPairs * 2, 0);
366 }
367 m_comm->GetColumnComm()->AllReduce(data, LibUtilities::ReduceSum);
368 }
369 m_recvMobInfo.clear();
370 for (int t = 0; t < totPairs; ++t)
371 {
372 int from, to;
373 unPack2Short(from, to, data[2 * t]);
374 if (to == m_rank)
375 {
376 m_recvMobInfo[(int)from] = data[2 * t + 1];
377 }
378 }
379}
void unPack2Short(int &a, int &b, const int &c)
void Pack2Short(const int &a, const int &b, int &c)

References m_columnRank, m_comm, m_mobilePts, m_rank, m_recvMobInfo, m_recvStatInfo, m_rowRank, CellMLToNektar.cellml_metadata::p, Pack2Short(), Nektar::LibUtilities::ReduceMax, Nektar::LibUtilities::ReduceSum, and unPack2Short().

Referenced by Nektar::SolverUtils::FilterLagrangianPoints::v_Initialise().

◆ SyncColumnComm()

void Nektar::SolverUtils::EvaluatePoints::SyncColumnComm ( )
protected

Definition at line 621 of file EvaluatePoints.cpp.

622{
624 {
625 return;
626 }
627 if (m_columnRank == 0)
628 {
629 for (const auto &p : m_mobilePts)
630 {
631 p.second->m_local[2] = p.second->m_global[2];
632 }
633 }
634 else
635 {
636 m_mobilePts.clear();
637 }
638 if (m_comm->GetColumnComm()->GetSize() == 1)
639 {
640 return;
641 }
642 int nPts = m_mobilePts.size();
643 m_comm->GetColumnComm()->AllReduce(nPts, LibUtilities::ReduceMax);
644 int nPackageSize = 1 + m_spacedim;
645 Array<OneD, NekDouble> data(nPts * nPackageSize, 0.), tmp;
646 if (m_columnRank == 0)
647 {
648 int count = 0;
649 for (const auto &p : m_mobilePts)
650 {
651 Pack2Int(p.first, p.second->m_eId, data[count++]);
652 Vmath::Vcopy(3, p.second->m_local, 1, tmp = data + count, 1);
653 count += 3;
654 }
655 }
656 m_comm->GetColumnComm()->AllReduce(data, LibUtilities::ReduceSum);
657 if (m_columnRank != 0)
658 {
659 int offset = 0;
660 for (int p = 0; p < nPts; ++p)
661 {
662 int pId, eId;
663 Array<OneD, NekDouble> globalCoords(m_spacedim, 0.);
664 unPack2Int(pId, eId, data[offset++]);
666 m_spacedim, 0, globalCoords);
667 m_mobilePts[pId]->SetLocalCoords(eId, data + offset);
668 offset += m_spacedim;
669 }
670 }
671}

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), m_columnRank, m_comm, m_isHomogeneous1D, m_mobilePts, m_spacedim, CellMLToNektar.cellml_metadata::p, Pack2Int(), Nektar::LibUtilities::ReduceMax, Nektar::LibUtilities::ReduceSum, unPack2Int(), and Vmath::Vcopy().

Referenced by PartitionMobilePts().

◆ unPack2Int()

void Nektar::SolverUtils::EvaluatePoints::unPack2Int ( int &  a,
int &  b,
const double &  d 
)
protected

Definition at line 62 of file EvaluatePoints.cpp.

63{
64 int *p = (int *)&d;
65 a = p[0];
66 b = p[1];
67}

References Nektar::UnitTests::d(), and CellMLToNektar.cellml_metadata::p.

Referenced by PartitionExchangeNonlocalPoints(), and SyncColumnComm().

◆ unPack2Short()

void Nektar::SolverUtils::EvaluatePoints::unPack2Short ( int &  a,
int &  b,
const int &  c 
)
protected

Definition at line 95 of file EvaluatePoints.cpp.

96{
97 short int *p = (short int *)&c;
98 a = (int)p[0];
99 b = (int)p[1];
100}

References CellMLToNektar.cellml_metadata::p.

Referenced by SetUpCommInfo().

◆ unPack3Int()

void Nektar::SolverUtils::EvaluatePoints::unPack3Int ( int &  a,
int &  b,
int &  c,
const double &  d 
)
protected

Definition at line 79 of file EvaluatePoints.cpp.

80{
81 int *p = (int *)&d;
82 unsigned short int *q = (unsigned short int *)&d;
83 a = p[0];
84 b = q[2];
85 c = q[3];
86}

References Nektar::UnitTests::d(), CellMLToNektar.cellml_metadata::p, and Nektar::UnitTests::q().

Referenced by PartitionExchangeNonlocalPoints(), PassMobilePhysToStatic(), and UpdateMobileCoords().

◆ UpdateMobileCoords()

void Nektar::SolverUtils::EvaluatePoints::UpdateMobileCoords ( std::map< int, Array< OneD, NekDouble > > &  globalCoords)
protected

Definition at line 688 of file EvaluatePoints.cpp.

690{
691 if (m_columnRank == 0 || !m_isHomogeneous1D)
692 {
693 for (auto &iter : globalCoords)
694 {
695 if (iter.first != m_rank)
696 {
697 m_comm->Send(iter.first, iter.second);
698 }
699 }
700 for (const auto &iter : m_recvStatInfo)
701 {
702 Array<OneD, NekDouble> data;
703 if (iter.first != m_rank)
704 {
705 data = Array<OneD, NekDouble>(iter.second * (m_spacedim + 1));
706 m_comm->Recv(iter.first, data);
707 }
708 else
709 {
710 data = globalCoords[iter.first];
711 }
712 for (int offset = 0; offset < data.size(); offset += m_spacedim)
713 {
714 int pid, mobt, stat;
715 unPack3Int(pid, mobt, stat, data[offset++]);
716 m_mobilePts[pid]->SetGlobalCoords(data + offset);
717 }
718 }
719 }
720}

References m_columnRank, m_comm, m_isHomogeneous1D, m_mobilePts, m_rank, m_recvStatInfo, m_spacedim, and unPack3Int().

Referenced by PassStaticCoordsToMobile().

◆ v_ModifyVelocity()

void Nektar::SolverUtils::EvaluatePoints::v_ModifyVelocity ( Array< OneD, NekDouble gcoords,
NekDouble  time,
Array< OneD, NekDouble vel 
)
protectedvirtual

Reimplemented in Nektar::SolverUtils::FilterLagrangianPoints.

Definition at line 266 of file EvaluatePoints.cpp.

270{
271}

Referenced by EvaluateMobilePhys().

Member Data Documentation

◆ m_columnRank

int Nektar::SolverUtils::EvaluatePoints::m_columnRank
protected

◆ m_comm

LibUtilities::CommSharedPtr Nektar::SolverUtils::EvaluatePoints::m_comm
protected

◆ m_defaultValueFunction

std::map<int, LibUtilities::EquationSharedPtr> Nektar::SolverUtils::EvaluatePoints::m_defaultValueFunction
protected

Definition at line 261 of file EvaluatePoints.h.

Referenced by GetDefaultValue(), and Initialise().

◆ m_isHomogeneous1D

bool Nektar::SolverUtils::EvaluatePoints::m_isHomogeneous1D
protected

◆ m_mobilePts

std::map<int, MobilePointSharedPtr> Nektar::SolverUtils::EvaluatePoints::m_mobilePts
protected

◆ m_nPhysics

int Nektar::SolverUtils::EvaluatePoints::m_nPhysics
protected

◆ m_rank

int Nektar::SolverUtils::EvaluatePoints::m_rank
protected

◆ m_recvMobInfo

std::map<int, int> Nektar::SolverUtils::EvaluatePoints::m_recvMobInfo
protected

◆ m_recvStatInfo

std::map<int, int> Nektar::SolverUtils::EvaluatePoints::m_recvStatInfo
protected

◆ m_rowRank

int Nektar::SolverUtils::EvaluatePoints::m_rowRank
protected

Definition at line 254 of file EvaluatePoints.h.

Referenced by Initialise(), PartitionExchangeNonlocalPoints(), and SetUpCommInfo().

◆ m_spacedim

int Nektar::SolverUtils::EvaluatePoints::m_spacedim
protected

◆ m_staticPts

StationaryPointsSharedPtr Nektar::SolverUtils::EvaluatePoints::m_staticPts
protected