Nektar++
Public Member Functions | Private Attributes | List of all members
Nektar::MultiRegions::InterfaceTrace Class Reference

#include <InterfaceMapDG.h>

Public Member Functions

 InterfaceTrace (const ExpListSharedPtr &trace, const SpatialDomains::InterfaceShPtr &interfaceShPtr)
 Constructor. More...
 
virtual ~InterfaceTrace ()=default
 Default destructor. More...
 
void SetCheckLocal (bool flag)
 
std::vector< Array< OneD, NekDouble > > GetMissingCoords ()
 Returns the missing coordinates vector. More...
 
SpatialDomains::InterfaceShPtr GetInterface ()
 Returns the interface object. More...
 
void CalcLocalMissing ()
 Calculates what coordinates on the interface are missing locally. More...
 
void FillLocalBwdTrace (Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
 Fills the Bwd trace by interpolating from the Fwd for local interfaces. More...
 
void FillRankBwdTrace (Array< OneD, NekDouble > &trace, Array< OneD, NekDouble > &Bwd)
 Fills the Bwd trace from partitioned trace. More...
 

Private Attributes

ExpListSharedPtr m_trace
 Trace expansion list. More...
 
SpatialDomains::InterfaceShPtr m_interface
 Local interface object. More...
 
bool m_checkLocal = false
 Flag whether the opposite side of the interface is present locally. More...
 
std::vector< Array< OneD, NekDouble > > m_missingCoords
 Vector of coordinates on interface missing from the other side locally. More...
 
std::map< int, std::pair< int, Array< OneD, NekDouble > > > m_foundLocalCoords
 Map of found coordinates present locally. More...
 
std::vector< int > m_mapMissingCoordToTrace
 Vector of indices corresponding to m_missingCoord locations in trace. More...
 

Detailed Description

Object for each interface present between two ranks, which are held in the InterfaceExchange object.

Definition at line 51 of file InterfaceMapDG.h.

Constructor & Destructor Documentation

◆ InterfaceTrace()

Nektar::MultiRegions::InterfaceTrace::InterfaceTrace ( const ExpListSharedPtr trace,
const SpatialDomains::InterfaceShPtr interfaceShPtr 
)

Constructor.

Definition at line 42 of file InterfaceMapDG.cpp.

45 : m_trace(trace), m_interface(interfaceShPtr)
46{
47}
ExpListSharedPtr m_trace
Trace expansion list.
SpatialDomains::InterfaceShPtr m_interface
Local interface object.

◆ ~InterfaceTrace()

virtual Nektar::MultiRegions::InterfaceTrace::~InterfaceTrace ( )
virtualdefault

Default destructor.

Member Function Documentation

◆ CalcLocalMissing()

void Nektar::MultiRegions::InterfaceTrace::CalcLocalMissing ( )

Calculates what coordinates on the interface are missing locally.

Calculates what coordinates on the interface are missing; i.e. aren't located in an edge from the other side of the interface on this partition. These are stored in the vector m_missingCoords, and another vector of the same index layout contains the location of that coordinate in the trace i.e. the 'offset + i' value. It checks first if the interface is flagged local, which denotes whether the opposite side is also present on this rank, if not then all coordinates are missing. If local then the same structure as CalcRankDistances is used to minimise computational cost.

Definition at line 245 of file InterfaceMapDG.cpp.

246{
247 auto childEdge = m_interface->GetEdge();
248 // If not flagged 'check local' then all points are missing
249 if (!m_checkLocal)
250 {
251 int cnt = 0;
252 for (auto &childId : childEdge)
253 {
254 auto childElmt = m_trace->GetExpFromGeomId(childId.first);
255 size_t nq = childElmt->GetTotPoints();
256 Array<OneD, NekDouble> xc(nq, 0.0), yc(nq, 0.0), zc(nq, 0.0);
257 childElmt->GetCoords(xc, yc, zc);
258 int offset =
259 m_trace->GetPhys_Offset(m_trace->GetElmtToExpId(childId.first));
260
261 for (int i = 0; i < nq; ++i, ++cnt)
262 {
263 Array<OneD, NekDouble> xs(3);
264 xs[0] = xc[i];
265 xs[1] = yc[i];
266 xs[2] = zc[i];
267
268 m_missingCoords.emplace_back(xs);
269 m_mapMissingCoordToTrace.emplace_back(offset + i);
270 }
271 }
272 }
273 // Otherwise we need to check to see what points the other side of the
274 // interface on this rank contains. This is done by looping over each
275 // quadrature point in all the edges on this side of the interface, then
276 // check if this quadrature point is also in any edges on the other side
277 // of the interface. First use MinMaxCheck as a cheap method to eliminate
278 // far away edges. If it passes MinMaxCheck perform the FindDistance search
279 // which numerically searches for the closest local coordinate on the other
280 // interface side edge and returns the cartesian distance from the search
281 // coordinate to the found local coordinate. If determined to be close
282 // enough then add to found local coordinates.
283 else
284 {
285 for (auto &childId : childEdge)
286 {
287 auto childElmt = m_trace->GetExpFromGeomId(childId.first);
288 size_t nq = childElmt->GetTotPoints();
289 Array<OneD, NekDouble> xc(nq, 0.0), yc(nq, 0.0), zc(nq, 0.0);
290 childElmt->GetCoords(xc, yc, zc);
291 int offset =
292 m_trace->GetPhys_Offset(m_trace->GetElmtToExpId(childId.first));
293
294 for (int i = 0; i < nq; ++i)
295 {
296 bool found = false;
297 Array<OneD, NekDouble> foundLocCoord;
298 Array<OneD, NekDouble> xs(3);
299 xs[0] = xc[i];
300 xs[1] = yc[i];
301 xs[2] = zc[i];
302
303 auto parentEdge = m_interface->GetOppInterface()->GetEdge();
304 for (auto &edge : parentEdge)
305 {
306 // First check if inside the edge bounding box
307 if (!edge.second->MinMaxCheck(xs))
308 {
309 continue;
310 }
311
312 NekDouble dist =
313 edge.second->FindDistance(xs, foundLocCoord);
315 {
316 found = true;
317 m_foundLocalCoords[offset + i] = std::make_pair(
318 edge.second->GetGlobalID(), foundLocCoord);
319 break;
320 }
321 }
322
323 if (!found)
324 {
325 m_missingCoords.emplace_back(xs);
326 m_mapMissingCoordToTrace.emplace_back(offset + i);
327 }
328 }
329 }
330 }
331
332 // If running in serial there shouldn't be any missing coordinates.
333 if (m_trace->GetComm()->IsSerial())
334 {
336 "Missing " + std::to_string(m_missingCoords.size()) +
337 " coordinates on interface ID " +
338 std::to_string(m_interface->GetId()) +
339 " linked to interface ID " +
340 std::to_string(m_interface->GetOppInterface()->GetId()) +
341 ". Check both sides of the interface line up.");
342 }
343}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
std::map< int, std::pair< int, Array< OneD, NekDouble > > > m_foundLocalCoords
Map of found coordinates present locally.
std::vector< int > m_mapMissingCoordToTrace
Vector of indices corresponding to m_missingCoord locations in trace.
std::vector< Array< OneD, NekDouble > > m_missingCoords
Vector of coordinates on interface missing from the other side locally.
bool m_checkLocal
Flag whether the opposite side of the interface is present locally.
static const NekDouble kFindDistanceMin
double NekDouble

References ASSERTL0, Nektar::NekConstants::kFindDistanceMin, m_checkLocal, m_foundLocalCoords, m_interface, m_mapMissingCoordToTrace, m_missingCoords, and m_trace.

◆ FillLocalBwdTrace()

void Nektar::MultiRegions::InterfaceTrace::FillLocalBwdTrace ( Array< OneD, NekDouble > &  Fwd,
Array< OneD, NekDouble > &  Bwd 
)

Fills the Bwd trace by interpolating from the Fwd for local interfaces.

Definition at line 472 of file InterfaceMapDG.cpp.

474{
475 // If flagged then fill trace from local coords
476 if (m_checkLocal)
477 {
478 for (auto &foundLocCoord : m_foundLocalCoords)
479 {
480 int traceId = m_trace->GetElmtToExpId(foundLocCoord.second.first);
481 Array<OneD, NekDouble> locCoord = foundLocCoord.second.second;
482
483 Array<OneD, NekDouble> edgePhys =
484 Fwd + m_trace->GetPhys_Offset(traceId);
485
486 Bwd[foundLocCoord.first] =
487 m_trace->GetExp(traceId)->StdPhysEvaluate(locCoord, edgePhys);
488 }
489 }
490}

References m_checkLocal, m_foundLocalCoords, and m_trace.

◆ FillRankBwdTrace()

void Nektar::MultiRegions::InterfaceTrace::FillRankBwdTrace ( Array< OneD, NekDouble > &  trace,
Array< OneD, NekDouble > &  Bwd 
)

Fills the Bwd trace from partitioned trace.

Definition at line 507 of file InterfaceMapDG.cpp.

509{
510 for (int i = 0; i < m_mapMissingCoordToTrace.size(); ++i)
511 {
512 if (!std::isnan(trace[i]))
513 {
514 Bwd[m_mapMissingCoordToTrace[i]] = trace[i];
515 }
516 }
517}

References m_mapMissingCoordToTrace.

◆ GetInterface()

SpatialDomains::InterfaceShPtr Nektar::MultiRegions::InterfaceTrace::GetInterface ( )
inline

Returns the interface object.

Definition at line 74 of file InterfaceMapDG.h.

75 {
76 return m_interface;
77 }

References m_interface.

◆ GetMissingCoords()

std::vector< Array< OneD, NekDouble > > Nektar::MultiRegions::InterfaceTrace::GetMissingCoords ( )
inline

Returns the missing coordinates vector.

Definition at line 68 of file InterfaceMapDG.h.

69 {
70 return m_missingCoords;
71 }

References m_missingCoords.

◆ SetCheckLocal()

void Nektar::MultiRegions::InterfaceTrace::SetCheckLocal ( bool  flag)
inline

Definition at line 62 of file InterfaceMapDG.h.

63 {
64 m_checkLocal = flag;
65 }

References m_checkLocal.

Member Data Documentation

◆ m_checkLocal

bool Nektar::MultiRegions::InterfaceTrace::m_checkLocal = false
private

Flag whether the opposite side of the interface is present locally.

Definition at line 94 of file InterfaceMapDG.h.

Referenced by CalcLocalMissing(), FillLocalBwdTrace(), and SetCheckLocal().

◆ m_foundLocalCoords

std::map<int, std::pair<int, Array<OneD, NekDouble> > > Nektar::MultiRegions::InterfaceTrace::m_foundLocalCoords
private

Map of found coordinates present locally.

Definition at line 98 of file InterfaceMapDG.h.

Referenced by CalcLocalMissing(), and FillLocalBwdTrace().

◆ m_interface

SpatialDomains::InterfaceShPtr Nektar::MultiRegions::InterfaceTrace::m_interface
private

Local interface object.

Definition at line 92 of file InterfaceMapDG.h.

Referenced by CalcLocalMissing(), and GetInterface().

◆ m_mapMissingCoordToTrace

std::vector<int> Nektar::MultiRegions::InterfaceTrace::m_mapMissingCoordToTrace
private

Vector of indices corresponding to m_missingCoord locations in trace.

Definition at line 100 of file InterfaceMapDG.h.

Referenced by CalcLocalMissing(), and FillRankBwdTrace().

◆ m_missingCoords

std::vector<Array<OneD, NekDouble> > Nektar::MultiRegions::InterfaceTrace::m_missingCoords
private

Vector of coordinates on interface missing from the other side locally.

Definition at line 96 of file InterfaceMapDG.h.

Referenced by CalcLocalMissing(), and GetMissingCoords().

◆ m_trace

ExpListSharedPtr Nektar::MultiRegions::InterfaceTrace::m_trace
private

Trace expansion list.

Definition at line 90 of file InterfaceMapDG.h.

Referenced by CalcLocalMissing(), and FillLocalBwdTrace().