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
37
39{
40
42 const ExpListSharedPtr &trace,
43 const SpatialDomains::InterfaceShPtr &interfaceShPtr)
44 : m_trace(trace), m_interface(interfaceShPtr)
45{
46}
47
48// Creates the objects that communicates trace values across the interface in
49// serial and parallel. Does this by checking all interfaces to see if they
50// are contained locally (on the current rank), and checks all other ranks for
51// the same interface ID and if contains the opposite interface side i.e.
52// 'left' or 'right'
55 const ExpListSharedPtr &trace)
56 : m_graph(meshGraph), m_movement(meshGraph->GetMovement()), m_trace(trace)
57{
58 auto comm = m_trace->GetComm()->GetSpaceComm();
59 auto interfaceCollection = m_movement->GetInterfaces();
60
61 // myIndxLR contains the info about what interface edges are present on
62 // current rank with each interface no, i, consisting of:
63 // [i] = indx
64 // [i + 1] = 0 (non), = 1 (left only), = 2 (right only), = 3 (both)
65 std::map<int, int> myIndxLRMap;
66
67 // localInterfaces contains a map of interface ID to a pair of interface
68 // traces, this pair is 'left' and 'right' interface and is used to
69 // construct and store the traces for interfaces present on the current rank
70 std::map<int, std::pair<InterfaceTraceSharedPtr, InterfaceTraceSharedPtr>>
71 localInterfaces;
72
73 // Map of m_localInterfaces vector to interface ID
74 Array<OneD, int> indxToInterfaceID(interfaceCollection.size());
75
76 // Loops over all interfaces and check if either side, 'left' or 'right'
77 // are empty on the current rank, if not empty then populate the local
78 // interface data structures and keep track using myIndxLRMap
79 size_t cnt = 0;
80 for (const auto &interface : interfaceCollection)
81 {
82 indxToInterfaceID[cnt] = interface.first.first;
83 myIndxLRMap[interface.first.first] = 0;
84
85 if (!interface.second->GetLeftInterface()->IsEmpty())
86 {
87 myIndxLRMap[interface.first.first] += 1;
88
89 localInterfaces[interface.first.first].first =
91 trace, interface.second->GetLeftInterface());
92 m_localInterfaces.emplace_back(
93 localInterfaces[interface.first.first].first);
94 }
95
96 if (!interface.second->GetRightInterface()->IsEmpty())
97 {
98 myIndxLRMap[interface.first.first] += 2;
99
100 localInterfaces[interface.first.first].second =
102 trace, interface.second->GetRightInterface());
103 m_localInterfaces.emplace_back(
104 localInterfaces[interface.first.first].second);
105 }
106
107 cnt++;
108 }
109
110 // Send num of interfaces size so all partitions can prepare buffers
111 int nRanks = comm->GetSize();
112
113 // Send all interface edges present to all partitions
114 Array<OneD, int> interfaceEdges(myIndxLRMap.size());
115 cnt = 0;
116 for (auto pres : myIndxLRMap)
117 {
118 interfaceEdges[cnt++] = pres.second;
119 }
120 Array<OneD, int> rankLocalInterfaceIds(myIndxLRMap.size() * nRanks, 0);
121 comm->AllGather(interfaceEdges, rankLocalInterfaceIds);
122
123 // Map of rank to vector of interface traces
124 std::map<int, std::vector<InterfaceTraceSharedPtr>> oppRankSharedInterface;
125
126 // Find what interface Ids match with other ranks, then check if opposite
127 // edge
128 size_t myRank = comm->GetRank();
129 size_t numInterfaces = interfaceCollection.size();
130 for (int i = 0; i < nRanks; ++i)
131 {
132 for (size_t j = 0; j < numInterfaces; ++j)
133 {
134 int otherId = indxToInterfaceID[j];
135
136 // otherCode represents for a specific interface ID what sides are
137 // present on rank i, 0=non, 1=left, 2=right, 3=both
138 int otherCode = rankLocalInterfaceIds[i * numInterfaces + j];
139
140 // myCode represents for a specific interface ID what sides are
141 // present on this rank/process 0=non, 1=left, 2=right, 3=both
142 int myCode = myIndxLRMap[otherId];
143
144 // Special case if checking current rank (this process)
145 if (i == myRank)
146 {
147 // If contains both edges locally then set check local to true
148 if (myCode == 3)
149 {
150 localInterfaces[otherId].first->SetCheckLocal(true);
151 localInterfaces[otherId].second->SetCheckLocal(true);
152 }
153
154 continue;
155 }
156
157 // Checks if this ranks 'left' matches any 'right' on other rank
158 if ((myCode == 1 && otherCode == 2) ||
159 (myCode == 1 && otherCode == 3) ||
160 (myCode == 3 && otherCode == 2))
161 {
162 oppRankSharedInterface[i].emplace_back(
163 localInterfaces[otherId].first);
164 }
165 // Checks if this ranks 'right' matches any 'left' on other rank
166 else if ((myCode == 2 && otherCode == 1) ||
167 (myCode == 2 && otherCode == 3) ||
168 (myCode == 3 && otherCode == 1))
169 {
170 oppRankSharedInterface[i].emplace_back(
171 localInterfaces[otherId].second);
172 }
173 // Checks if this ranks 'both' matches any 'both' on other rank
174 else if (myCode == 3 && otherCode == 3)
175 {
176 oppRankSharedInterface[i].emplace_back(
177 localInterfaces[otherId].first);
178 oppRankSharedInterface[i].emplace_back(
179 localInterfaces[otherId].second);
180 }
181 }
182 }
183
184 // Create individual interface exchange objects (each object is rank ->
185 // rank) and contains a vector of interfaceTrace objects
186 for (auto &rank : oppRankSharedInterface)
187 {
188 m_exchange.emplace_back(
190 m_movement, m_trace, comm, rank));
191 }
192
193 // Find missing coordinates on interface from other side
195}
196
198{
200 timer.Start();
201
202 auto comm = m_trace->GetComm();
203 auto zones = m_movement->GetZones();
204
205 LibUtilities::Timer timer2;
206 timer2.Start();
207 for (auto &interfaceTrace : m_localInterfaces)
208 {
209 interfaceTrace->CalcLocalMissing(m_movement);
210 }
211 timer2.Stop();
212 timer2.AccumulateRegion("InterfaceMapDG::ExchangeCoords local", 1);
213
214 // If parallel communication is needed
215 if (!m_exchange.empty())
216 {
217 LibUtilities::Timer timer3;
218 timer3.Start();
219
220 auto requestSend = comm->CreateRequest(m_exchange.size());
221 auto requestRecv = comm->CreateRequest(m_exchange.size());
222
223 for (int i = 0; i < m_exchange.size(); ++i)
224 {
225 m_exchange[i]->RankFillSizes(requestSend, requestRecv, i);
226 }
227 comm->WaitAll(requestSend);
228 comm->WaitAll(requestRecv);
229
230 for (int i = 0; i < m_exchange.size(); ++i)
231 {
232 m_exchange[i]->SendMissing(requestSend, requestRecv, i);
233 }
234 comm->WaitAll(requestSend);
235 comm->WaitAll(requestRecv);
236
237 for (auto &i : m_exchange)
238 {
239 i->CalcRankDistances();
240 }
241
242 timer3.Stop();
243 timer3.AccumulateRegion("InterfaceMapDG::ExchangeCoords parallel", 1);
244 }
245
246 m_movement->GetCoordExchangeFlag() = false;
247
248 timer.Stop();
249 timer.AccumulateRegion("InterfaceMapDG::ExchangeCoords");
250}
251
252/**
253 * Calculates what coordinates on the interface are missing; i.e. aren't located
254 * in an edge from the other side of the interface on this partition. These are
255 * stored in the vector #m_missingCoords, and another vector of the same index
256 * layout contains the location of that coordinate in the trace i.e. the
257 * 'offset + i' value. It checks first if the interface is flagged local, which
258 * denotes whether the opposite side is also present on this rank, if not then
259 * all coordinates are missing. If local then the same structure as
260 * CalcRankDistances is used to minimise computational cost.
261 */
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 {
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;
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}
399
400// This fills m_sendSize and m_recvSize with the correct sizes for the
401// communication lengths for each rank when communicating the traces. I.e.
402// however many quadrature points on the opposite side of the interface to
403// send, and however many on this side to receive.
406 LibUtilities::CommRequestSharedPtr &requestRecv, int requestNum)
407{
408
409 // Recalc size is the number of interfaces missing points that need to be
410 // communicated.
411 int recalcSize = 0;
412 for (auto &localInterface : m_interfaceTraces)
413 {
414 recalcSize +=
415 m_zones[localInterface->GetInterface()->GetId()]->GetMoved();
416 }
417
418 m_sendSize[m_rank] = Array<OneD, int>(recalcSize);
419 m_recvSize[m_rank] = Array<OneD, int>(recalcSize);
420
421 int cnt = 0;
422 for (auto &interface : m_interfaceTraces)
423 {
424 if (m_zones[interface->GetInterface()->GetId()]->GetMoved())
425 {
426 // Calculates the number of missing coordinates to be communicated
427 // This size is 3 times the number of missing coordinates to allow
428 // for 3D points
429 m_sendSize[m_rank][cnt++] =
430 interface->GetMissingCoords().size() * 3;
431 }
432 }
433
434 // Uses non-blocking send/recvs to send communication sizes to other ranks
436 requestSend, requestNum);
438 requestRecv, requestNum);
439}
440
441// Uses the calculated send/recv sizes to send all missing coordinates to
442// other ranks.
445 LibUtilities::CommRequestSharedPtr &requestRecv, int requestNum)
446{
447 m_totSendSize[m_rank] = std::accumulate(m_sendSize[m_rank].begin(),
448 m_sendSize[m_rank].end(), 0);
449 m_totRecvSize[m_rank] = std::accumulate(m_recvSize[m_rank].begin(),
450 m_recvSize[m_rank].end(), 0);
451
454
455 int cnt = 0;
456 for (auto &interface : m_interfaceTraces)
457 {
458 if (m_zones[interface->GetInterface()->GetId()]->GetMoved())
459 {
460 auto missing = interface->GetMissingCoords();
461 for (auto coord : missing)
462 {
463 // Points are organised in blocks of three
464 // x coord, then y coord, then z coord.
465 for (int k = 0; k < 3; ++k, ++cnt)
466 {
467 m_send[cnt] = coord[k];
468 }
469 }
470 }
471 }
472
473 // Uses non-blocking send/recvs to send missing points to other ranks
474 m_comm->Isend(m_rank, m_send, m_totSendSize[m_rank], requestSend,
475 requestNum);
476 m_comm->Irecv(m_rank, m_recv, m_totRecvSize[m_rank], requestRecv,
477 requestNum);
478}
479
480/**
481 * Performs the trace exchange across interfaces
482 * 1) Calculate and send the Fwd trace to other ranks with pairwise comms
483 * 2) While waiting for the send/recv to complete we fill the local interfaces
484 * Bwd trace
485 * 3) Fill the remaining Bwd trace with the received trace data from other ranks
486 */
489{
490 if (m_movement->GetCoordExchangeFlag())
491 {
493 }
494 auto comm = m_trace->GetComm();
495 // If no parallel exchange needed we only fill the local traces
496 if (m_exchange.empty())
497 {
498
499 // Fill local interface traces
500 for (auto &m_localInterface : m_localInterfaces)
501 {
502 m_localInterface->FillLocalBwdTrace(Fwd, Bwd);
503 }
504 }
505 else
506 {
507 auto requestSend = comm->CreateRequest(m_exchange.size());
508 auto requestRecv = comm->CreateRequest(m_exchange.size());
509 for (int i = 0; i < m_exchange.size(); ++i)
510 {
511 m_exchange[i]->SendFwdTrace(requestSend, requestRecv, i, Fwd);
512 }
513
514 // Fill local interface traces
515 for (auto &m_localInterface : m_localInterfaces)
516 {
517 m_localInterface->FillLocalBwdTrace(Fwd, Bwd);
518 }
519
520 comm->WaitAll(requestSend);
521 comm->WaitAll(requestRecv);
522
523 // Fill communicated interface traces
524 for (auto &i : m_exchange)
525 {
526 i->FillRankBwdTraceExchange(Bwd);
527 }
528 }
529}
530
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}
550
552{
553 int cnt = 0;
554 for (int i = 0; i < m_interfaceTraces.size(); ++i)
555 {
556 Array<OneD, NekDouble> traceTmp(m_sendSize[m_rank][i] / 3, 0.0);
557 for (int j = 0; j < m_sendSize[m_rank][i] / 3; ++j, ++cnt)
558 {
559 traceTmp[j] = m_recvTrace[cnt];
560 }
561
562 m_interfaceTraces[i]->FillRankBwdTrace(traceTmp, Bwd);
563 }
564}
565
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}
577
580 LibUtilities::CommRequestSharedPtr &requestRecv, int requestNum,
582{
584 Array<OneD, NekDouble>(m_totSendSize[m_rank] / 3, std::nan(""));
586 Array<OneD, NekDouble>(m_totRecvSize[m_rank] / 3, std::nan(""));
587
588 for (auto &i : m_foundRankCoords[m_rank])
589 {
590 int traceId = m_trace->GetElmtToExpId(i.second.first.second);
591 Array<OneD, NekDouble> locCoord = i.second.second;
592
593 Array<OneD, NekDouble> edgePhys =
594 Fwd + m_trace->GetPhys_Offset(traceId);
595
596 m_sendTrace[i.first] =
597 m_trace->GetExp(traceId)->StdPhysEvaluate(locCoord, edgePhys);
598 }
599
600 m_comm->Isend(m_rank, m_sendTrace, m_sendTrace.size(), requestSend,
601 requestNum);
602 m_comm->Irecv(m_rank, m_recvTrace, m_recvTrace.size(), requestRecv,
603 requestNum);
604}
605
606/**
607 * Check coords in m_recv from other rank to see if present on this rank, and
608 * then populates m_foundRankCoords if found.
609 * - Loops over all interface traces present in the exchange, i.e. every
610 * interface needed on this rank-to-rank communication object
611 * - Then loops over all the quadrature points missing from the other side in
612 * that interface using global coordinates
613 * - Loop over each local edge and use MinMaxCheck on the point as a fast
614 * comparison on whether the point lies on the edge
615 * - If true then use the more expensive FindDistance function to find the
616 * closest local coordinate of the global 'missing' point on the local edge
617 * and the distance from the edge that corresponds to
618 * - If distance is less than 5e-5 save found coordinate in m_foundRankCoords
619 */
620
622{
623 // Clear old found coordinates
624 auto foundRankCoordsCopy = m_foundRankCoords[m_rank];
626 .clear(); // @TODO: This may cause problems with 2 interfaces and one is
627 // fixed? With the skip below.
628
629 Array<OneD, int> disp(m_recvSize[m_rank].size() + 1, 0.0);
630 std::partial_sum(m_recvSize[m_rank].begin(), m_recvSize[m_rank].end(),
631 &disp[1]);
632
633 for (int i = 0; i < m_interfaceTraces.size(); ++i)
634 {
635 if (!m_zones[m_interfaceTraces[i]->GetInterface()->GetId()]->GetMoved())
636 {
637 // If zone is not moved then skip
638 continue;
639 }
640
641 auto localEdge = m_interfaceTraces[i]->GetInterface()->GetEdge();
642
643 for (int j = disp[i]; j < disp[i + 1]; j += 3)
644 {
645 Array<OneD, NekDouble> foundLocCoord;
647 xs[0] = m_recv[j];
648 xs[1] = m_recv[j + 1];
649 xs[2] = m_recv[j + 2];
650
651 // First search the edge the point was found in last timestep
652 if ((foundRankCoordsCopy.find(j / 3) !=
653 foundRankCoordsCopy.end()) &&
654 (foundRankCoordsCopy[j / 3].first.first == i))
655 {
656 auto edge = m_interfaceTraces[i]->GetInterface()->GetEdge(
657 foundRankCoordsCopy[j / 3].first.second);
658 NekDouble dist = edge->FindDistance(xs, foundLocCoord);
659
661 {
662 m_foundRankCoords[m_rank][j / 3] = std::make_pair(
663 std::make_pair(i, edge->GetGlobalID()), foundLocCoord);
664 continue;
665 }
666 }
667
668 for (auto &edge : localEdge)
669 {
670 // First check if inside the edge bounding box
671 if (!edge.second->MinMaxCheck(xs))
672 {
673 continue;
674 }
675
676 NekDouble dist = edge.second->FindDistance(xs, foundLocCoord);
677
679 {
680 m_foundRankCoords[m_rank][j / 3] = std::make_pair(
681 std::make_pair(i, edge.second->GetGlobalID()),
682 foundLocCoord);
683 break;
684 }
685 }
686 }
687 }
688}
689
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}
744
745} // namespace Nektar::MultiRegions
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
void AccumulateRegion(std::string, int iolevel=0)
Accumulate elapsed time for a region.
Definition: Timer.cpp:70
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, 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.
std::map< int, std::map< int, std::pair< std::pair< int, int >, Array< OneD, NekDouble > > > > m_foundRankCoords
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.
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.
void CalcLocalMissing(SpatialDomains::MovementSharedPtr movement)
Calculates what coordinates on the interface are missing locally.
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 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< Movement > MovementSharedPtr
Definition: MeshGraph.h:177
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:174
double NekDouble