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

#include <InterfaceMapDG.h>

Public Member Functions

 ~InterfaceMapDG ()=default
 Default destructor. More...
 
 InterfaceMapDG (const SpatialDomains::MeshGraphSharedPtr &graph, const ExpListSharedPtr &trace)
 
void ExchangeTrace (Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
 Perform the trace exchange between processors, given the forwards and backwards spaces. More...
 
void ExchangeCoords ()
 Perform the coordinate exchange between processors. This is where the missing coordinates on the interface are found and sent to all other processors on the other side of that interface so the matching ranks can be found. More...
 

Private Attributes

SpatialDomains::MeshGraphSharedPtr m_graph
 Mesh associated with this expansion list. More...
 
SpatialDomains::MovementSharedPtr m_movement
 Movement object associated with the non-conformal interfaces. More...
 
std::vector< InterfaceTraceSharedPtrm_localInterfaces
 Interface sides present on current process. More...
 
const ExpListSharedPtr m_trace
 Trace expansion list. More...
 
std::vector< InterfaceExchangeSharedPtrm_exchange
 Vector of interface exchanges, i.e. every rank-to-rank comm needed. More...
 

Detailed Description

Implements the communication patterns to allow for exchange of information across non-conformal interfaces and across different partitions. Holds all the InterfaceExchange objects in m_exchange.

Definition at line 221 of file InterfaceMapDG.h.

Constructor & Destructor Documentation

◆ ~InterfaceMapDG()

Nektar::MultiRegions::InterfaceMapDG::~InterfaceMapDG ( )
default

Default destructor.

◆ InterfaceMapDG()

Nektar::MultiRegions::InterfaceMapDG::InterfaceMapDG ( const SpatialDomains::MeshGraphSharedPtr graph,
const ExpListSharedPtr trace 
)

Sets up the InterfaceExchange objects stored in m_exchange, each object is rank -> rank and contains a vector of InterfaceTrace objects corresponding to shared interfaces between those ranks.

Definition at line 52 of file InterfaceMapDG.cpp.

55 : m_graph(meshGraph), m_movement(meshGraph->GetMovement()), m_trace(trace)
56{
57 auto comm = m_trace->GetComm()->GetSpaceComm();
58 auto interfaceCollection = m_movement->GetInterfaces();
59
60 // myIndxLR contains the info about what interface edges are present on
61 // current rank with each interface no, i, consisting of:
62 // [i] = indx
63 // [i + 1] = 0 (non), = 1 (left only), = 2 (right only), = 3 (both)
64 std::map<int, int> myIndxLRMap;
65
66 // localInterfaces contains a map of interface ID to a pair of interface
67 // traces, this pair is 'left' and 'right' interface and is used to
68 // construct and store the traces for interfaces present on the current rank
69 std::map<int, std::pair<InterfaceTraceSharedPtr, InterfaceTraceSharedPtr>>
70 localInterfaces;
71
72 // Map of m_localInterfaces vector to interface ID
73 Array<OneD, int> indxToInterfaceID(interfaceCollection.size());
74
75 // Loops over all interfaces and check if either side, 'left' or 'right'
76 // are empty on the current rank, if not empty then populate the local
77 // interface data structures and keep track using myIndxLRMap
78 size_t cnt = 0;
79 for (const auto &interface : interfaceCollection)
80 {
81 indxToInterfaceID[cnt] = interface.first.first;
82 myIndxLRMap[interface.first.first] = 0;
83
84 if (!interface.second->GetLeftInterface()->IsEmpty())
85 {
86 myIndxLRMap[interface.first.first] += 1;
87
88 localInterfaces[interface.first.first].first =
90 trace, interface.second->GetLeftInterface());
91 m_localInterfaces.emplace_back(
92 localInterfaces[interface.first.first].first);
93 }
94
95 if (!interface.second->GetRightInterface()->IsEmpty())
96 {
97 myIndxLRMap[interface.first.first] += 2;
98
99 localInterfaces[interface.first.first].second =
101 trace, interface.second->GetRightInterface());
102 m_localInterfaces.emplace_back(
103 localInterfaces[interface.first.first].second);
104 }
105
106 cnt++;
107 }
108
109 // Send num of interfaces size so all partitions can prepare buffers
110 int nRanks = comm->GetSize();
111
112 // Send all interface edges present to all partitions
113 Array<OneD, int> interfaceEdges(myIndxLRMap.size());
114 cnt = 0;
115 for (auto pres : myIndxLRMap)
116 {
117 interfaceEdges[cnt++] = pres.second;
118 }
119 Array<OneD, int> rankLocalInterfaceIds(myIndxLRMap.size() * nRanks, 0);
120 comm->AllGather(interfaceEdges, rankLocalInterfaceIds);
121
122 // Map of rank to vector of interface traces
123 std::map<int, std::vector<InterfaceTraceSharedPtr>> oppRankSharedInterface;
124
125 // Find what interface Ids match with other ranks, then check if opposite
126 // edge
127 size_t myRank = comm->GetRank();
128 size_t numInterfaces = interfaceCollection.size();
129 for (int i = 0; i < nRanks; ++i)
130 {
131 for (size_t j = 0; j < numInterfaces; ++j)
132 {
133 int otherId = indxToInterfaceID[j];
134
135 // otherCode represents for a specific interface ID what sides are
136 // present on rank i, 0=non, 1=left, 2=right, 3=both
137 int otherCode = rankLocalInterfaceIds[i * numInterfaces + j];
138
139 // myCode represents for a specific interface ID what sides are
140 // present on this rank/process 0=non, 1=left, 2=right, 3=both
141 int myCode = myIndxLRMap[otherId];
142
143 // Special case if checking current rank (this process)
144 if (i == myRank)
145 {
146 // If contains both edges locally then set check local to true
147 if (myCode == 3)
148 {
149 localInterfaces[otherId].first->SetCheckLocal(true);
150 localInterfaces[otherId].second->SetCheckLocal(true);
151 }
152
153 continue;
154 }
155
156 // Checks if this ranks 'left' matches any 'right' on other rank
157 if ((myCode == 1 && otherCode == 2) ||
158 (myCode == 1 && otherCode == 3) ||
159 (myCode == 3 && otherCode == 2))
160 {
161 oppRankSharedInterface[i].emplace_back(
162 localInterfaces[otherId].first);
163 }
164 // Checks if this ranks 'right' matches any 'left' on other rank
165 else if ((myCode == 2 && otherCode == 1) ||
166 (myCode == 2 && otherCode == 3) ||
167 (myCode == 3 && otherCode == 1))
168 {
169 oppRankSharedInterface[i].emplace_back(
170 localInterfaces[otherId].second);
171 }
172 // Checks if this ranks 'both' matches any 'both' on other rank
173 else if (myCode == 3 && otherCode == 3)
174 {
175 oppRankSharedInterface[i].emplace_back(
176 localInterfaces[otherId].first);
177 oppRankSharedInterface[i].emplace_back(
178 localInterfaces[otherId].second);
179 }
180 }
181 }
182
183 // Create individual interface exchange objects (each object is rank ->
184 // rank) and contains a vector of interfaceTrace objects
185 for (auto &rank : oppRankSharedInterface)
186 {
187 m_exchange.emplace_back(
189 m_movement, m_trace, comm, rank));
190 }
191
192 // Find missing coordinates on interface from other side
194}
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
SpatialDomains::MeshGraphSharedPtr m_graph
Mesh associated with this expansion list.
const ExpListSharedPtr m_trace
Trace expansion list.
std::vector< InterfaceExchangeSharedPtr > m_exchange
Vector of interface exchanges, i.e. every rank-to-rank comm needed.
void ExchangeCoords()
Perform the coordinate exchange between processors. This is where the missing coordinates on the inte...
SpatialDomains::MovementSharedPtr m_movement
Movement object associated with the non-conformal interfaces.
std::vector< InterfaceTraceSharedPtr > m_localInterfaces
Interface sides present on current process.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), ExchangeCoords(), m_exchange, m_localInterfaces, m_movement, and m_trace.

Member Function Documentation

◆ ExchangeCoords()

void Nektar::MultiRegions::InterfaceMapDG::ExchangeCoords ( )

Perform the coordinate exchange between processors. This is where the missing coordinates on the interface are found and sent to all other processors on the other side of that interface so the matching ranks can be found.

Definition at line 196 of file InterfaceMapDG.cpp.

197{
198 auto comm = m_trace->GetComm()->GetSpaceComm();
199 auto zones = m_movement->GetZones();
200
201 for (auto &interfaceTrace : m_localInterfaces)
202 {
203 interfaceTrace->CalcLocalMissing();
204 }
205
206 // If parallel communication is needed
207 if (!m_exchange.empty())
208 {
209 auto requestSend = comm->CreateRequest(m_exchange.size());
210 auto requestRecv = comm->CreateRequest(m_exchange.size());
211
212 for (int i = 0; i < m_exchange.size(); ++i)
213 {
214 m_exchange[i]->RankFillSizes(requestSend, requestRecv, i);
215 }
216 comm->WaitAll(requestSend);
217 comm->WaitAll(requestRecv);
218
219 for (int i = 0; i < m_exchange.size(); ++i)
220 {
221 m_exchange[i]->SendMissing(requestSend, requestRecv, i);
222 }
223 comm->WaitAll(requestSend);
224 comm->WaitAll(requestRecv);
225
226 for (auto &i : m_exchange)
227 {
228 i->CalcRankDistances();
229 }
230 }
231}

References m_exchange, m_localInterfaces, m_movement, and m_trace.

Referenced by InterfaceMapDG().

◆ ExchangeTrace()

void Nektar::MultiRegions::InterfaceMapDG::ExchangeTrace ( Array< OneD, NekDouble > &  Fwd,
Array< OneD, NekDouble > &  Bwd 
)

Perform the trace exchange between processors, given the forwards and backwards spaces.

Parameters
FwdLocal forwards space of the trace (which will be sent)
BwdLocal backwards space of the trace (which will receive contributions)

Performs the trace exchange across interfaces 1) Calculate and send the Fwd trace to other ranks with pairwise comms 2) While waiting for the send/recv to complete we fill the local interfaces Bwd trace 3) Fill the remaining Bwd trace with the received trace data from other ranks

Definition at line 430 of file InterfaceMapDG.cpp.

432{
433 auto comm = m_trace->GetComm()->GetSpaceComm();
434
435 // If no parallel exchange needed we only fill the local traces
436 if (m_exchange.empty())
437 {
438 // Fill local interface traces
439 for (auto &m_localInterface : m_localInterfaces)
440 {
441 m_localInterface->FillLocalBwdTrace(Fwd, Bwd);
442 }
443 }
444 else
445 {
446 auto requestSend = comm->CreateRequest(m_exchange.size());
447 auto requestRecv = comm->CreateRequest(m_exchange.size());
448 for (int i = 0; i < m_exchange.size(); ++i)
449 {
450 m_exchange[i]->SendFwdTrace(requestSend, requestRecv, i, Fwd);
451 }
452
453 // Fill local interface traces
454 for (auto &m_localInterface : m_localInterfaces)
455 {
456 m_localInterface->FillLocalBwdTrace(Fwd, Bwd);
457 }
458
459 comm->WaitAll(requestSend);
460 comm->WaitAll(requestRecv);
461
462 // Fill communicated interface traces
463 for (auto &i : m_exchange)
464 {
465 i->FillRankBwdTraceExchange(Bwd);
466 }
467 }
468}

References m_exchange, m_localInterfaces, and m_trace.

Member Data Documentation

◆ m_exchange

std::vector<InterfaceExchangeSharedPtr> Nektar::MultiRegions::InterfaceMapDG::m_exchange
private

Vector of interface exchanges, i.e. every rank-to-rank comm needed.

Definition at line 264 of file InterfaceMapDG.h.

Referenced by ExchangeCoords(), ExchangeTrace(), and InterfaceMapDG().

◆ m_graph

SpatialDomains::MeshGraphSharedPtr Nektar::MultiRegions::InterfaceMapDG::m_graph
private

Mesh associated with this expansion list.

Definition at line 256 of file InterfaceMapDG.h.

◆ m_localInterfaces

std::vector<InterfaceTraceSharedPtr> Nektar::MultiRegions::InterfaceMapDG::m_localInterfaces
private

Interface sides present on current process.

Definition at line 260 of file InterfaceMapDG.h.

Referenced by ExchangeCoords(), ExchangeTrace(), and InterfaceMapDG().

◆ m_movement

SpatialDomains::MovementSharedPtr Nektar::MultiRegions::InterfaceMapDG::m_movement
private

Movement object associated with the non-conformal interfaces.

Definition at line 258 of file InterfaceMapDG.h.

Referenced by ExchangeCoords(), and InterfaceMapDG().

◆ m_trace

const ExpListSharedPtr Nektar::MultiRegions::InterfaceMapDG::m_trace
private

Trace expansion list.

Definition at line 262 of file InterfaceMapDG.h.

Referenced by ExchangeCoords(), ExchangeTrace(), and InterfaceMapDG().