Nektar++
InterfaceMapDG.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File InterfaceMapDG.cpp
4//
5// For more information, please see: http://www.nektar.info
6//
7// The MIT License
8//
9// Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10// Department of Aeronautics, Imperial College London (UK), and Scientific
11// Computing and Imaging Institute, University of Utah (USA).
12//
13// Permission is hereby granted, free of charge, to any person obtaining a
14// copy of this software and associated documentation files (the "Software"),
15// to deal in the Software without restriction, including without limitation
16// the rights to use, copy, modify, merge, publish, distribute, sublicense,
17// and/or sell copies of the Software, and to permit persons to whom the
18// Software is furnished to do so, subject to the following conditions:
19//
20// The above copyright notice and this permission notice shall be included
21// in all copies or substantial portions of the Software.
22//
23// THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24// OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25// FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26// THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27// LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28// FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29// DEALINGS IN THE SOFTWARE.
30//
31// Description: MPI communication for interfaces
32//
33///////////////////////////////////////////////////////////////////////////////
34
35#include "InterfaceMapDG.h"
36
38{
39
41 const ExpListSharedPtr &trace,
42 const SpatialDomains::InterfaceShPtr &interfaceShPtr)
43 : m_trace(trace), m_interface(interfaceShPtr)
44{
45}
46
47// Creates the objects that communicates trace values across the interface in
48// serial and parallel. Does this by checking all interfaces to see if they
49// are contained locally (on the current rank), and checks all other ranks for
50// the same interface ID and if contains the opposite interface side i.e.
51// 'left' or 'right'
54 const ExpListSharedPtr &trace)
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}
195
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}
232
233/**
234 * Calculates what coordinates on the interface are missing; i.e. aren't located
235 * in an edge from the other side of the interface on this partition. These are
236 * stored in the vector #m_missingCoords, and another vector of the same index
237 * layout contains the location of that coordinate in the trace i.e. the
238 * 'offset + i' value. It checks first if the interface is flagged local, which
239 * denotes whether the opposite side is also present on this rank, if not then
240 * all coordinates are missing. If local then the same structure as
241 * CalcRankDistances is used to minimise computational cost.
242 */
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 {
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;
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}
342
343// This fills m_sendSize and m_recvSize with the correct sizes for the
344// communication lengths for each rank when communicating the traces. I.e.
345// however many quadrature points on the opposite side of the interface to
346// send, and however many on this side to receive.
349 LibUtilities::CommRequestSharedPtr &requestRecv, int requestNum)
350{
351
352 // Recalc size is the number of interfaces missing points that need to be
353 // communicated.
354 int recalcSize = 0;
355 for (auto &localInterface : m_interfaceTraces)
356 {
357 recalcSize +=
358 m_zones[localInterface->GetInterface()->GetId()]->GetMoved();
359 }
360
361 m_sendSize[m_rank] = Array<OneD, int>(recalcSize);
362 m_recvSize[m_rank] = Array<OneD, int>(recalcSize);
363
364 int cnt = 0;
365 for (auto &interface : m_interfaceTraces)
366 {
367 if (m_zones[interface->GetInterface()->GetId()]->GetMoved())
368 {
369 // Calculates the number of missing coordinates to be communicated
370 // This size is 3 times the number of missing coordinates to allow
371 // for 3D points
372 m_sendSize[m_rank][cnt++] =
373 interface->GetMissingCoords().size() * 3;
374 }
375 }
376
377 // Uses non-blocking send/recvs to send communication sizes to other ranks
379 requestSend, requestNum);
381 requestRecv, requestNum);
382}
383
384// Uses the calculated send/recv sizes to send all missing coordinates to
385// other ranks.
388 LibUtilities::CommRequestSharedPtr &requestRecv, int requestNum)
389{
390 m_totSendSize[m_rank] = std::accumulate(m_sendSize[m_rank].begin(),
391 m_sendSize[m_rank].end(), 0);
392 m_totRecvSize[m_rank] = std::accumulate(m_recvSize[m_rank].begin(),
393 m_recvSize[m_rank].end(), 0);
394
397
398 int cnt = 0;
399 for (auto &interface : m_interfaceTraces)
400 {
401 if (m_zones[interface->GetInterface()->GetId()]->GetMoved())
402 {
403 auto missing = interface->GetMissingCoords();
404 for (auto coord : missing)
405 {
406 // Points are organised in blocks of three
407 // x coord, then y coord, then z coord.
408 for (int k = 0; k < 3; ++k, ++cnt)
409 {
410 m_send[cnt] = coord[k];
411 }
412 }
413 }
414 }
415
416 // Uses non-blocking send/recvs to send missing points to other ranks
417 m_comm->Isend(m_rank, m_send, m_totSendSize[m_rank], requestSend,
418 requestNum);
419 m_comm->Irecv(m_rank, m_recv, m_totRecvSize[m_rank], requestRecv,
420 requestNum);
421}
422
423/**
424 * Performs the trace exchange across interfaces
425 * 1) Calculate and send the Fwd trace to other ranks with pairwise comms
426 * 2) While waiting for the send/recv to complete we fill the local interfaces
427 * Bwd trace
428 * 3) Fill the remaining Bwd trace with the received trace data from other ranks
429 */
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}
469
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}
489
491{
492 int cnt = 0;
493 for (int i = 0; i < m_interfaceTraces.size(); ++i)
494 {
495 Array<OneD, NekDouble> traceTmp(m_sendSize[m_rank][i] / 3, 0.0);
496 for (int j = 0; j < m_sendSize[m_rank][i] / 3; ++j, ++cnt)
497 {
498 traceTmp[j] = m_recvTrace[cnt];
499 }
500
501 m_interfaceTraces[i]->FillRankBwdTrace(traceTmp, Bwd);
502 }
503}
504
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}
516
519 LibUtilities::CommRequestSharedPtr &requestRecv, int requestNum,
521{
523 Array<OneD, NekDouble>(m_totSendSize[m_rank] / 3, std::nan(""));
525 Array<OneD, NekDouble>(m_totRecvSize[m_rank] / 3, std::nan(""));
526
527 for (auto &i : m_foundRankCoords[m_rank])
528 {
529 int traceId = m_trace->GetElmtToExpId(i.second.first);
530 Array<OneD, NekDouble> locCoord = i.second.second;
531
532 Array<OneD, NekDouble> edgePhys =
533 Fwd + m_trace->GetPhys_Offset(traceId);
534
535 m_sendTrace[i.first] =
536 m_trace->GetExp(traceId)->StdPhysEvaluate(locCoord, edgePhys);
537 }
538
539 m_comm->Isend(m_rank, m_sendTrace, m_sendTrace.size(), requestSend,
540 requestNum);
541 m_comm->Irecv(m_rank, m_recvTrace, m_recvTrace.size(), requestRecv,
542 requestNum);
543}
544
545/**
546 * Check coords in m_recv from other rank to see if present on this rank, and
547 * then populates m_foundRankCoords if found.
548 * - Loops over all interface traces present in the exchange, i.e. every
549 * interface needed on this rank-to-rank communication object
550 * - Then loops over all the quadrature points missing from the other side in
551 * that interface using global coordinates
552 * - Loop over each local edge and use MinMaxCheck on the point as a fast
553 * comparison on whether the point lies on the edge
554 * - If true then use the more expensive FindDistance function to find the
555 * closest local coordinate of the global 'missing' point on the local edge
556 * and the distance from the edge that corresponds to
557 * - If distance is less than 5e-5 save found coordinate in m_foundRankCoords
558 */
559
561{
562 Array<OneD, int> disp(m_recvSize[m_rank].size() + 1, 0.0);
563 std::partial_sum(m_recvSize[m_rank].begin(), m_recvSize[m_rank].end(),
564 &disp[1]);
565
566 for (int i = 0; i < m_interfaceTraces.size(); ++i)
567 {
568 auto localEdge = m_interfaceTraces[i]->GetInterface()->GetEdge();
569
570 for (int j = disp[i]; j < disp[i + 1]; j += 3)
571 {
572 Array<OneD, NekDouble> foundLocCoord;
574 xs[0] = m_recv[j];
575 xs[1] = m_recv[j + 1];
576 xs[2] = m_recv[j + 2];
577
578 for (auto &edge : localEdge)
579 {
580 // First check if inside the edge bounding box
581 if (!edge.second->MinMaxCheck(xs))
582 {
583 continue;
584 }
585
586 NekDouble dist = edge.second->FindDistance(xs, foundLocCoord);
587
589 {
590 m_foundRankCoords[m_rank][j / 3] = std::make_pair(
591 edge.second->GetGlobalID(), foundLocCoord);
592 break;
593 }
594 }
595 }
596 }
597}
598
599} // namespace Nektar::MultiRegions
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
General purpose memory allocation routines with the ability to allocate from thread specific memory p...
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
std::map< int, Array< OneD, int > > m_sendSize
Map of rank to array of size of send buffer for each interface.
void RankFillSizes(LibUtilities::CommRequestSharedPtr &requestSend, LibUtilities::CommRequestSharedPtr &requestRecv, int requestNum)
std::map< int, SpatialDomains::ZoneBaseShPtr > m_zones
Map of zone IDs to zone bases.
const ExpListSharedPtr m_trace
Trace expansion list.
void FillRankBwdTraceExchange(Array< OneD, NekDouble > &Bwd)
const std::vector< InterfaceTraceSharedPtr > m_interfaceTraces
Vector of interface traces i.e. every interface side present on m_rank.
const LibUtilities::CommSharedPtr m_comm
Communicator.
void SendFwdTrace(LibUtilities::CommRequestSharedPtr &requestSend, LibUtilities::CommRequestSharedPtr &requestRecv, int requestNum, Array< OneD, NekDouble > &Fwd)
Array< OneD, NekDouble > m_send
Send buffer for coord exchange.
std::map< int, std::map< int, std::pair< int, Array< OneD, NekDouble > > > > m_foundRankCoords
std::map< int, int > m_totSendSize
Map of rank to total size of send buffer for all interfaces.
void CalcRankDistances()
Populates m_foundRankCoords using the FindDistance function.
std::map< int, int > m_totRecvSize
Map of rank to total size of receive buffer for all interfaces.
std::map< int, Array< OneD, int > > m_recvSize
Map of rank to array of size of receive buffer for each interface.
Array< OneD, NekDouble > m_recv
Receive buffer for coord exchange.
void SendMissing(LibUtilities::CommRequestSharedPtr &requestSend, LibUtilities::CommRequestSharedPtr &requestRecv, int requestNum)
Array< OneD, NekDouble > m_recvTrace
Receive buffer for trace exchange.
Array< OneD, NekDouble > m_sendTrace
Send buffer for trace exchange.
const ExpListSharedPtr m_trace
Trace expansion list.
void ExchangeTrace(Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
Perform the trace exchange between processors, given the forwards and backwards spaces.
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.
InterfaceMapDG(const SpatialDomains::MeshGraphSharedPtr &graph, const ExpListSharedPtr &trace)
std::vector< InterfaceTraceSharedPtr > m_localInterfaces
Interface sides present on current process.
ExpListSharedPtr m_trace
Trace expansion list.
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.
SpatialDomains::InterfaceShPtr m_interface
Local interface object.
std::vector< Array< OneD, NekDouble > > m_missingCoords
Vector of coordinates on interface missing from the other side locally.
InterfaceTrace(const ExpListSharedPtr &trace, const SpatialDomains::InterfaceShPtr &interfaceShPtr)
Constructor.
void FillRankBwdTrace(Array< OneD, NekDouble > &trace, Array< OneD, NekDouble > &Bwd)
Fills the Bwd trace from partitioned trace.
void CalcLocalMissing()
Calculates what coordinates on the interface are missing locally.
void FillLocalBwdTrace(Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
Fills the Bwd trace by interpolating from the Fwd for local interfaces.
bool m_checkLocal
Flag whether the opposite side of the interface is present locally.
std::shared_ptr< CommRequest > CommRequestSharedPtr
Definition: Comm.h:84
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
static const NekDouble kFindDistanceMin
std::shared_ptr< Interface > InterfaceShPtr
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:174
double NekDouble