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

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

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

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

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

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

◆ m_localInterfaces

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

Interface sides present on current process.

Definition at line 262 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 260 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 264 of file InterfaceMapDG.h.

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