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 49 of file InterfaceMapDG.h.

Constructor & Destructor Documentation

◆ InterfaceTrace()

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

Constructor.

Definition at line 40 of file InterfaceMapDG.cpp.

43 : m_trace(trace), m_interface(interfaceShPtr)
44{
45}
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 243 of file InterfaceMapDG.cpp.

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

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

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

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

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 92 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 96 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 90 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 98 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 94 of file InterfaceMapDG.h.

Referenced by CalcLocalMissing(), and GetMissingCoords().

◆ m_trace

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

Trace expansion list.

Definition at line 88 of file InterfaceMapDG.h.

Referenced by CalcLocalMissing(), and FillLocalBwdTrace().