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 (SpatialDomains::MovementSharedPtr movement)
 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...
 
void DomainCheck (Array< OneD, NekDouble > &gloCoord, SpatialDomains::MovementSharedPtr movement)
 Check whether the coordiniates in the original domain. 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 49 of file InterfaceMapDG.h.

Constructor & Destructor Documentation

◆ InterfaceTrace()

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

Constructor.

Definition at line 41 of file InterfaceMapDG.cpp.

44 : m_trace(trace), m_interface(interfaceShPtr)
45{
46}
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 ( SpatialDomains::MovementSharedPtr  movement)

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 262 of file InterfaceMapDG.cpp.

264{
265 // Nuke old missing/found
266 m_missingCoords.clear();
268
269 // Store copy of found coords to optimise search before clearing
270 auto foundLocalCoordsCopy = m_foundLocalCoords;
271 m_foundLocalCoords.clear();
272
273 auto childEdge = m_interface->GetEdge();
274 // If not flagged 'check local' then all points are missing
275 if (!m_checkLocal)
276 {
277 int cnt = 0;
278 for (auto &childId : childEdge)
279 {
280 auto childElmt = m_trace->GetExpFromGeomId(childId.first);
281 size_t nq = childElmt->GetTotPoints();
282 Array<OneD, NekDouble> xc(nq, 0.0), yc(nq, 0.0), zc(nq, 0.0);
283 childElmt->GetCoords(xc, yc, zc);
284 int offset =
285 m_trace->GetPhys_Offset(m_trace->GetElmtToExpId(childId.first));
286
287 for (int i = 0; i < nq; ++i, ++cnt)
288 {
289 Array<OneD, NekDouble> xs(3);
290 xs[0] = xc[i];
291 xs[1] = yc[i];
292 xs[2] = zc[i];
293
294 if (movement->GetMovedFlag() && movement->GetTranslateFlag())
295 {
296 DomainCheck(xs, movement);
297 }
298
299 m_missingCoords.emplace_back(xs);
300 m_mapMissingCoordToTrace.emplace_back(offset + i);
301 }
302 }
303 }
304 // Otherwise we need to check to see what points the other side of the
305 // interface on this rank contains. This is done by looping over each
306 // quadrature point in all the edges on this side of the interface, then
307 // check if this quadrature point is also in any edges on the other side
308 // of the interface. First use MinMaxCheck as a cheap method to eliminate
309 // far away edges. If it passes MinMaxCheck perform the FindDistance search
310 // which numerically searches for the closest local coordinate on the other
311 // interface side edge and returns the cartesian distance from the search
312 // coordinate to the found local coordinate. If determined to be close
313 // enough then add to found local coordinates.
314 else
315 {
316 for (auto &childId : childEdge)
317 {
318 auto childElmt = m_trace->GetExpFromGeomId(childId.first);
319 size_t nq = childElmt->GetTotPoints();
320 Array<OneD, NekDouble> xc(nq, 0.0), yc(nq, 0.0), zc(nq, 0.0);
321 childElmt->GetCoords(xc, yc, zc);
322 int offset =
323 m_trace->GetPhys_Offset(m_trace->GetElmtToExpId(childId.first));
324
325 for (int i = 0; i < nq; ++i)
326 {
327 bool found = false;
328 Array<OneD, NekDouble> foundLocCoord;
329 Array<OneD, NekDouble> xs(3);
330 xs[0] = xc[i];
331 xs[1] = yc[i];
332 xs[2] = zc[i];
333
334 if (movement->GetMovedFlag() && movement->GetTranslateFlag())
335 {
336 DomainCheck(xs, movement);
337 }
338
339 // First search the edge the point was found in last timestep
340 if (foundLocalCoordsCopy.find(offset + i) !=
341 foundLocalCoordsCopy.end())
342 {
343 auto edge = m_interface->GetOppInterface()->GetEdge(
344 foundLocalCoordsCopy[offset + i].first);
345 NekDouble dist = edge->FindDistance(xs, foundLocCoord);
346
348 {
349 m_foundLocalCoords[offset + i] =
350 std::make_pair(edge->GetGlobalID(), foundLocCoord);
351 continue;
352 }
353 }
354
355 // If not in last timestep edge then loop over all interface
356 // edges
357 auto parentEdge = m_interface->GetOppInterface()->GetEdge();
358 for (auto &edge : parentEdge)
359 {
360 // First check if inside the edge bounding box
361 if (!edge.second->MinMaxCheck(xs))
362 {
363 continue;
364 }
365
366 NekDouble dist =
367 edge.second->FindDistance(xs, foundLocCoord);
369 {
370 found = true;
371 m_foundLocalCoords[offset + i] = std::make_pair(
372 edge.second->GetGlobalID(), foundLocCoord);
373 break;
374 }
375 }
376
377 if (!found)
378 {
379 m_missingCoords.emplace_back(xs);
380 m_mapMissingCoordToTrace.emplace_back(offset + i);
381 }
382 }
383 }
384 }
385
386 // If running in serial there shouldn't be any missing coordinates.
387 // Unless we flag to ignore this check for this specific interface.
388 if (m_trace->GetComm()->IsSerial() && !m_interface->GetSkipCoordCheck())
389 {
391 "Missing " + std::to_string(m_missingCoords.size()) +
392 " coordinates on interface ID " +
393 std::to_string(m_interface->GetId()) +
394 " linked to interface ID " +
395 std::to_string(m_interface->GetOppInterface()->GetId()) +
396 ". Check both sides of the interface line up.");
397 }
398}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
void DomainCheck(Array< OneD, NekDouble > &gloCoord, SpatialDomains::MovementSharedPtr movement)
Check whether the coordiniates in the original domain.
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, DomainCheck(), Nektar::NekConstants::kFindDistanceMin, m_checkLocal, m_foundLocalCoords, m_interface, m_mapMissingCoordToTrace, m_missingCoords, and m_trace.

◆ DomainCheck()

void Nektar::MultiRegions::InterfaceTrace::DomainCheck ( Array< OneD, NekDouble > &  gloCoord,
SpatialDomains::MovementSharedPtr  movement 
)

Check whether the coordiniates in the original domain.

Definition at line 690 of file InterfaceMapDG.cpp.

692{
693 // Get current interface id and opposite interface id
694 int id = m_interface->GetId();
695 int oppid = m_interface->GetOppInterface()->GetId();
696
697 // Get displacement of the current interface and the opposite interface
698 auto zones = movement->GetZones();
699 auto disp = zones[id]->v_GetDisp();
700 auto oppdisp = zones[oppid]->v_GetDisp();
701
702 // Get domain length and domian box
703 auto box = movement->GetDomainBox();
704 auto length = movement->GetDomainLength();
705
706 NekDouble ExactMove = 0.0;
707 NekDouble ExactOppMove = 0.0;
708 int tmp = 0;
709
710 // update coordinates by add or minus domian length
711 for (int i = 0; i < disp.size(); ++i)
712 {
713 // map coordinate back to the original mesh
714 ExactMove = std::fmod(disp[i], length[i]);
715 ExactOppMove = std::fmod(oppdisp[i], length[i]);
716 tmp = oppdisp[i] / length[i];
717 gloCoord[i] = gloCoord[i] - disp[i] + ExactMove;
718
719 // update coordinates for perioidic interface
720 if (oppdisp[i] - disp[i] < 0)
721 {
722 if (gloCoord[i] > box[i + 3] + ExactOppMove)
723 {
724 gloCoord[i] = gloCoord[i] - length[i] + tmp * length[i];
725 }
726 else
727 {
728 gloCoord[i] = gloCoord[i] + tmp * length[i];
729 }
730 }
731 else
732 {
733 if (gloCoord[i] < box[i] + ExactOppMove)
734 {
735 gloCoord[i] = gloCoord[i] + length[i] + tmp * length[i];
736 }
737 else
738 {
739 gloCoord[i] = gloCoord[i] + tmp * length[i];
740 }
741 }
742 }
743}

References m_interface.

Referenced by CalcLocalMissing().

◆ 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 531 of file InterfaceMapDG.cpp.

533{
534 // If flagged then fill trace from local coords
535 if (m_checkLocal)
536 {
537 for (auto &foundLocCoord : m_foundLocalCoords)
538 {
539 int traceId = m_trace->GetElmtToExpId(foundLocCoord.second.first);
540 Array<OneD, NekDouble> locCoord = foundLocCoord.second.second;
541
542 Array<OneD, NekDouble> edgePhys =
543 Fwd + m_trace->GetPhys_Offset(traceId);
544
545 Bwd[foundLocCoord.first] =
546 m_trace->GetExp(traceId)->StdPhysEvaluate(locCoord, edgePhys);
547 }
548 }
549}

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 566 of file InterfaceMapDG.cpp.

568{
569 for (int i = 0; i < m_mapMissingCoordToTrace.size(); ++i)
570 {
571 if (!std::isnan(trace[i]))
572 {
573 Bwd[m_mapMissingCoordToTrace[i]] = trace[i];
574 }
575 }
576}

References m_mapMissingCoordToTrace.

◆ GetInterface()

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

Returns the interface object.

Definition at line 72 of file InterfaceMapDG.h.

73 {
74 return m_interface;
75 }

References m_interface.

◆ GetMissingCoords()

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

Returns the missing coordinates vector.

Definition at line 66 of file InterfaceMapDG.h.

67 {
68 return m_missingCoords;
69 }

References m_missingCoords.

◆ SetCheckLocal()

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

Definition at line 60 of file InterfaceMapDG.h.

61 {
62 m_checkLocal = flag;
63 }

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 95 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 99 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 93 of file InterfaceMapDG.h.

Referenced by CalcLocalMissing(), DomainCheck(), 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 101 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 97 of file InterfaceMapDG.h.

Referenced by CalcLocalMissing(), and GetMissingCoords().

◆ m_trace

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

Trace expansion list.

Definition at line 91 of file InterfaceMapDG.h.

Referenced by CalcLocalMissing(), and FillLocalBwdTrace().