Nektar++
Loading...
Searching...
No Matches
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, SpatialDomains::MovementSharedPtr movementShPtr)
 Constructor.
 
virtual ~InterfaceTrace ()=default
 Default destructor.
 
void SetCheckLocal (bool flag)
 
std::vector< Array< OneD, NekDouble > > GetMissingCoords ()
 Returns the missing coordinates vector.
 
std::map< int, std::pair< int, Array< OneD, NekDouble > > > GetLocalCoords ()
 
std::map< int, std::pair< int, NekDouble > > GetRotAngle ()
 
SpatialDomains::InterfaceShPtr GetInterface ()
 Returns the interface object.
 
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.
 
void RotLocalBwdTrace (Array< OneD, Array< OneD, NekDouble > > &Bwd, const int &dim)
 Rotates the interface velocity for sector rotation.
 
void RotLocalBwdDeriveTrace (TensorOfArray3D< NekDouble > &Bwd, const int &dim)
 Rotates the interface derivative velocity for sector rotation.
 
void DeriveRotate (TensorOfArray3D< NekDouble > &Bwd, const int dir, const NekDouble angle)
 Rotate the 3D vector field in Bwd around the axis 'dir' by 'angle'.
 
void FillRankBwdTrace (Array< OneD, NekDouble > &trace, Array< OneD, NekDouble > &Bwd)
 Fills the Bwd trace from partitioned trace.
 
void DomainCheck (Array< OneD, NekDouble > &gloCoord, NekDouble &angle)
 Update the coordiniates after movement.
 

Private Attributes

ExpListSharedPtr m_trace
 Trace expansion list.
 
SpatialDomains::InterfaceShPtr m_interface
 Local interface object.
 
SpatialDomains::MovementSharedPtr m_movement
 Movement object associated with the non-conformal interfaces.
 
bool m_checkLocal = false
 Flag whether the opposite side of the interface is present locally.
 
std::vector< Array< OneD, NekDouble > > m_missingCoords
 Vector of coordinates on interface missing from the other side locally.
 
std::map< int, std::pair< int, Array< OneD, NekDouble > > > m_foundLocalCoords
 Map of found coordinates present locally.
 
std::map< int, std::pair< int, NekDouble > > m_phyRotAngle
 Map of phy rotate angle for sector rotate interfaces.
 
std::vector< int > m_mapMissingCoordToTrace
 Vector of indices corresponding to m_missingCoord locations in trace.
 

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,
SpatialDomains::MovementSharedPtr  movementShPtr 
)

Constructor.

Definition at line 41 of file InterfaceMapDG.cpp.

45 : m_trace(trace), m_interface(interfaceShPtr), m_movement(movementShPtr)
46{
47}
ExpListSharedPtr m_trace
Trace expansion list.
SpatialDomains::InterfaceShPtr m_interface
Local interface object.
SpatialDomains::MovementSharedPtr m_movement
Movement object associated with the non-conformal interfaces.

◆ ~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 263 of file InterfaceMapDG.cpp.

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 // Zone id
273 int id = m_interface->GetId();
274 NekDouble angle = 0.0;
275
276 auto childEdge = m_interface->GetEdge();
277 // If not flagged 'check local' then all points are missing
278 if (!m_checkLocal)
279 {
280 int cnt = 0;
281 for (auto &childId : childEdge)
282 {
283 auto childElmt = m_trace->GetExpFromGeomId(childId.first);
284 size_t nq = childElmt->GetTotPoints();
285 Array<OneD, NekDouble> xc(nq, 0.0), yc(nq, 0.0), zc(nq, 0.0);
286 childElmt->GetCoords(xc, yc, zc);
287 int offset =
288 m_trace->GetPhys_Offset(m_trace->GetElmtToExpId(childId.first));
289 for (int i = 0; i < nq; ++i, ++cnt)
290 {
291 Array<OneD, NekDouble> xs(3);
292 xs[0] = xc[i];
293 xs[1] = yc[i];
294 xs[2] = zc[i];
295
296 // If mesh has moved, need to update corresponding points on
297 // interface
298 if (m_movement->GetMovedFlag())
299 {
300 DomainCheck(xs, angle);
301
302 // Record the rotation angle for sector rotation
303 if (m_movement->GetSectorRotateFlag())
304 {
305 m_phyRotAngle[offset + i] = std::make_pair(id, angle);
306 }
307 }
308
309 m_missingCoords.emplace_back(xs);
310 m_mapMissingCoordToTrace.emplace_back(offset + i);
311 }
312 }
313 }
314 // Otherwise we need to check to see what points the other side of the
315 // interface on this rank contains. This is done by looping over each
316 // quadrature point in all the edges on this side of the interface, then
317 // check if this quadrature point is also in any edges on the other side
318 // of the interface. First use MinMaxCheck as a cheap method to eliminate
319 // far away edges. If it passes MinMaxCheck perform the FindDistance search
320 // which numerically searches for the closest local coordinate on the other
321 // interface side edge and returns the cartesian distance from the search
322 // coordinate to the found local coordinate. If determined to be close
323 // enough then add to found local coordinates.
324 else
325 {
326 for (auto &childId : childEdge)
327 {
328 auto childElmt = m_trace->GetExpFromGeomId(childId.first);
329 size_t nq = childElmt->GetTotPoints();
330 Array<OneD, NekDouble> xc(nq, 0.0), yc(nq, 0.0), zc(nq, 0.0);
331 childElmt->GetCoords(xc, yc, zc);
332 int offset =
333 m_trace->GetPhys_Offset(m_trace->GetElmtToExpId(childId.first));
334 for (int i = 0; i < nq; ++i)
335 {
336 bool found = false;
337 Array<OneD, NekDouble> foundLocCoord;
338 Array<OneD, NekDouble> xs(3);
339 xs[0] = xc[i];
340 xs[1] = yc[i];
341 xs[2] = zc[i];
342
343 // If mesh has moved, need to update corresponding points on
344 // interface
345 if (m_movement->GetMovedFlag())
346 {
347 DomainCheck(xs, angle);
348
349 // Record the rotation angle for sector rotation
350 if (m_movement->GetSectorRotateFlag())
351 {
352 m_phyRotAngle[offset + i] = std::make_pair(id, angle);
353 }
354 }
355
356 // First search the edge the point was found in last timestep
357 if (foundLocalCoordsCopy.find(offset + i) !=
358 foundLocalCoordsCopy.end())
359 {
360 auto edge = m_interface->GetOppInterface()->GetEdge(
361 foundLocalCoordsCopy[offset + i].first);
362 NekDouble dist = edge->FindDistance(xs, foundLocCoord);
364 {
365 m_foundLocalCoords[offset + i] =
366 std::make_pair(edge->GetGlobalID(), foundLocCoord);
367 continue;
368 }
369 }
370
371 // If not in last timestep edge then loop over all interface
372 // edges
373 auto parentEdge = m_interface->GetOppInterface()->GetEdge();
374 for (auto &edge : parentEdge)
375 {
376 // First check if inside the edge bounding box
377 if (!edge.second->MinMaxCheck(xs))
378 {
379 continue;
380 }
381
382 NekDouble dist =
383 edge.second->FindDistance(xs, foundLocCoord);
385 {
386 found = true;
387 m_foundLocalCoords[offset + i] = std::make_pair(
388 edge.second->GetGlobalID(), foundLocCoord);
389 break;
390 }
391 }
392
393 if (!found)
394 {
395 m_missingCoords.emplace_back(xs);
396 m_mapMissingCoordToTrace.emplace_back(offset + i);
397 }
398 }
399 }
400 }
401
402 // If running in serial there shouldn't be any missing coordinates.
403 // Unless we flag to ignore this check for this specific interface.
404 if (m_trace->GetComm()->IsSerial() && !m_interface->GetSkipCoordCheck())
405 {
407 "Missing " + std::to_string(m_missingCoords.size()) +
408 " coordinates on interface ID " +
409 std::to_string(m_interface->GetId()) +
410 " linked to interface ID " +
411 std::to_string(m_interface->GetOppInterface()->GetId()) +
412 ". Check both sides of the interface line up.");
413 }
414}
#define ASSERTL0(condition, msg)
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.
void DomainCheck(Array< OneD, NekDouble > &gloCoord, NekDouble &angle)
Update the coordiniates after movement.
std::map< int, std::pair< int, NekDouble > > m_phyRotAngle
Map of phy rotate angle for sector rotate interfaces.
bool m_checkLocal
Flag whether the opposite side of the interface is present locally.
static const NekDouble kFindDistanceMin

References ASSERTL0, DomainCheck(), Nektar::NekConstants::kFindDistanceMin, m_checkLocal, m_foundLocalCoords, m_interface, m_mapMissingCoordToTrace, m_missingCoords, m_movement, m_phyRotAngle, and m_trace.

◆ DeriveRotate()

void Nektar::MultiRegions::InterfaceTrace::DeriveRotate ( TensorOfArray3D< NekDouble > &  Bwd,
const int  dir,
const NekDouble  angle 
)

Rotate the 3D vector field in Bwd around the axis 'dir' by 'angle'.

Definition at line 993 of file InterfaceMapDG.cpp.

995{
996 // Precompute trigonometric values using correct angle input
997 const NekDouble cosA = cos(angle);
998 const NekDouble sinA = sin(angle);
999
1000 // Define rotation matrix R (3x3) using Array<OneD, NekDouble> in
1001 // **column-major order**
1002 Array<OneD, NekDouble> R(9, 0.0);
1003
1004 if (dir == 0) // Rotation around X-axis
1005 {
1006 R[0] = 1.0;
1007 R[3] = 0.0;
1008 R[6] = 0.0;
1009 R[1] = 0.0;
1010 R[4] = cosA;
1011 R[7] = -sinA;
1012 R[2] = 0.0;
1013 R[5] = sinA;
1014 R[8] = cosA;
1015 }
1016 else if (dir == 1) // Rotation around Y-axis
1017 {
1018 R[0] = cosA;
1019 R[3] = 0.0;
1020 R[6] = sinA;
1021 R[1] = 0.0;
1022 R[4] = 1.0;
1023 R[7] = 0.0;
1024 R[2] = -sinA;
1025 R[5] = 0.0;
1026 R[8] = cosA;
1027 }
1028 else if (dir == 2) // Rotation around Z-axis
1029 {
1030 R[0] = cosA;
1031 R[3] = -sinA;
1032 R[6] = 0.0;
1033 R[1] = sinA;
1034 R[4] = cosA;
1035 R[7] = 0.0;
1036 R[2] = 0.0;
1037 R[5] = 0.0;
1038 R[8] = 1.0;
1039 }
1040 else
1041 {
1042 NEKERROR(
1044 "Invalid rotation axis for first-order derivative transformation.");
1045 }
1046
1047 // Allocate gradient tensor (3x3) in **column-major order**
1048 Array<OneD, NekDouble> gradU(9, 0.0);
1049
1050 // Load original velocity gradient tensor into **column-major** format
1051 gradU[0] = Bwd[0][1][0]; // ∂u_x/∂x
1052 gradU[3] = Bwd[1][1][0]; // ∂u_x/∂y
1053 gradU[6] = Bwd[2][1][0]; // ∂u_x/∂z
1054
1055 gradU[1] = Bwd[0][2][0]; // ∂u_y/∂x
1056 gradU[4] = Bwd[1][2][0]; // ∂u_y/∂y
1057 gradU[7] = Bwd[2][2][0]; // ∂u_y/∂z
1058
1059 gradU[2] = Bwd[0][3][0]; // ∂u_z/∂x
1060 gradU[5] = Bwd[1][3][0]; // ∂u_z/∂y
1061 gradU[8] = Bwd[2][3][0]; // ∂u_z/∂z
1062
1063 // Allocate intermediate matrix R * gradU in **column-major order**
1064 Array<OneD, NekDouble> R_gradU(9, 0.0);
1065
1066 // Perform matrix multiplication R * gradU using Blas::Dgemm in
1067 // **column-major order**
1068 Blas::Dgemm('N', 'N', 3, 3, 3, 1.0, R.data(), 3, gradU.data(), 3, 0.0,
1069 R_gradU.data(), 3);
1070
1071 // Allocate final rotated gradient (R * gradU) * R^T in **column-major
1072 // order**
1073 Array<OneD, NekDouble> gradU_rotated(9, 0.0);
1074
1075 // Perform matrix multiplication (R * gradU) * R^T using Blas::Dgemm in
1076 // **column-major order**
1077 Blas::Dgemm('N', 'T', 3, 3, 3, 1.0, R_gradU.data(), 3, R.data(), 3, 0.0,
1078 gradU_rotated.data(), 3);
1079
1080 // Store rotated gradients back in Bwd in **column-major order**
1081 Bwd[0][1][0] = gradU_rotated[0]; // ∂u_x'/∂x'
1082 Bwd[1][1][0] = gradU_rotated[3]; // ∂u_x'/∂y'
1083 Bwd[2][1][0] = gradU_rotated[6]; // ∂u_x'/∂z'
1084
1085 Bwd[0][2][0] = gradU_rotated[1]; // ∂u_y'/∂x'
1086 Bwd[1][2][0] = gradU_rotated[4]; // ∂u_y'/∂y'
1087 Bwd[2][2][0] = gradU_rotated[7]; // ∂u_y'/∂z'
1088
1089 Bwd[0][3][0] = gradU_rotated[2]; // ∂u_z'/∂x'
1090 Bwd[1][3][0] = gradU_rotated[5]; // ∂u_z'/∂y'
1091 Bwd[2][3][0] = gradU_rotated[8]; // ∂u_z'/∂z'
1092}
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
static void Dgemm(const char &transa, const char &transb, const int &m, const int &n, const int &k, const double &alpha, const double *a, const int &lda, const double *b, const int &ldb, const double &beta, double *c, const int &ldc)
BLAS level 3: Matrix-matrix multiply C = A x B where op(A)[m x k], op(B)[k x n], C[m x n] DGEMM perfo...
Definition Blas.hpp:383

References Blas::Dgemm(), Nektar::ErrorUtil::efatal, and NEKERROR.

Referenced by RotLocalBwdDeriveTrace().

◆ DomainCheck()

void Nektar::MultiRegions::InterfaceTrace::DomainCheck ( Array< OneD, NekDouble > &  gloCoord,
NekDouble angle 
)

Update the coordiniates after movement.

Definition at line 810 of file InterfaceMapDG.cpp.

812{
813 // Get current interface id and opposite interface id
814 int id = m_interface->GetId();
815 int oppid = m_interface->GetOppInterface()->GetId();
816
817 // Get displacement of the current interface and the opposite interface
818 auto zones = m_movement->GetZones();
819
820 if (m_movement->GetTranslateFlag())
821 {
822 auto disp = zones[id]->v_GetDisp();
823 auto oppdisp = zones[oppid]->v_GetDisp();
824 Array<OneD, NekDouble> box(6);
825 Array<OneD, NekDouble> length(3);
826
827 // Get zone length and domian box from the translate zone
828 if (zones[id]->GetMovementType() ==
830 {
831 auto zone = std::static_pointer_cast<SpatialDomains::ZoneTranslate>(
832 zones[id]);
833 box = zone->GetZoneBox();
834 length = zone->GetZoneLength();
835 }
836 else if (zones[oppid]->GetMovementType() ==
838 {
839 auto zone = std::static_pointer_cast<SpatialDomains::ZoneTranslate>(
840 zones[oppid]);
841 box = zone->GetZoneBox();
842 length = zone->GetZoneLength();
843 }
844
845 NekDouble ExactMove = 0.0;
846 NekDouble ExactOppMove = 0.0;
847 int tmp = 0;
848
849 // update coordinates by add or minus domian length
850 for (int i = 0; i < disp.size(); ++i)
851 {
852
853 // map coordinate back to the original mesh
854 ExactMove = std::fmod(disp[i], length[i]);
855 ExactOppMove = std::fmod(oppdisp[i], length[i]);
856 tmp = oppdisp[i] / length[i];
857 gloCoord[i] = gloCoord[i] - disp[i] + ExactMove;
858
859 // update coordinates for perioidic interface
860 if (oppdisp[i] - disp[i] < 0)
861 {
862 if (gloCoord[i] > box[i + 3] + ExactOppMove + 1e-8)
863 {
864 gloCoord[i] = gloCoord[i] - length[i] + tmp * length[i];
865 }
866 else
867 {
868 gloCoord[i] = gloCoord[i] + tmp * length[i];
869 }
870 }
871 else
872 {
873 if (gloCoord[i] < box[i] + ExactOppMove - 1e-8)
874 {
875 gloCoord[i] = gloCoord[i] + length[i] + tmp * length[i];
876 }
877 else
878 {
879 gloCoord[i] = gloCoord[i] + tmp * length[i];
880 }
881 }
882 }
883 }
884
885 if (m_movement->GetSectorRotateFlag())
886 {
887 // Using rotation zone information to update the coordinates
888 auto zone =
889 std::static_pointer_cast<SpatialDomains::ZoneRotate>(zones[id]);
890 auto zone2 =
891 std::static_pointer_cast<SpatialDomains::ZoneRotate>(zones[oppid]);
892
893 // Get the rotation angle of the current interface and the opposite
894 // interface
895 auto angle = zone->v_GetAngle();
896 auto oppAngle = zone2->v_GetAngle();
897
898 // Get the sector angle
899 NekDouble SectorAngle = zone->GetSectorAngle();
900
901 NekDouble ExactOppMove = 0.0;
902 RotateAngle = 0.0;
903
904 int tmp1 = 0;
905 int tmp2 = 0;
906 int tmp3 = 0;
907
908 // Get the base coordinate of the sector
909 Array<OneD, NekDouble> base = zone->GetSectorBase();
910
911 // Normalize angles to [0, 2*PI)
912 angle = std::fmod(angle, 2 * M_PI);
913 if (angle < 0)
914 {
915 angle += 2 * M_PI;
916 }
917 oppAngle = std::fmod(oppAngle, 2 * M_PI);
918 if (oppAngle < 0)
919 {
920 oppAngle += 2 * M_PI;
921 }
922
923 // Map coordinate back to the original mesh
924 ExactOppMove = std::fmod(oppAngle, SectorAngle);
925 tmp1 = static_cast<int>(std::floor(angle / SectorAngle));
926 tmp2 = static_cast<int>(std::floor(oppAngle / SectorAngle));
927
928 // Compute the angle between the current coordinate and the base
929 // coordinate
930 auto axis = zone->GetAxis();
931 NekDouble diff = 0.0;
932 if (axis[0] == 0 && axis[1] == 0 && axis[2] == 1)
933 {
934 diff = atan2(gloCoord[1], gloCoord[0]) - atan2(base[1], base[0]);
935 }
936 else if (axis[0] == 1 && axis[1] == 0 && axis[2] == 0)
937 {
938 diff = atan2(gloCoord[2], gloCoord[1]) - atan2(base[2], base[1]);
939 }
940 else if (axis[0] == 0 && axis[1] == 1 && axis[2] == 0)
941 {
942 diff = atan2(gloCoord[0], gloCoord[2]) - atan2(base[0], base[2]);
943 }
944 else
945 {
946 NEKERROR(ErrorUtil::efatal, "Invalid rotation axis.");
947 }
948
949 // Normalize angles to [0, 2*PI)
950 if (diff < 0)
951 {
952 diff += 2 * M_PI;
953 }
954 tmp3 = static_cast<int>(std::floor(diff / SectorAngle));
955
956 // Compute the exact difference after mapping the angle back to the
957 // original mesh
958 NekDouble Exactdiff = std::fmod(diff, SectorAngle);
959
960 // update coordinates for perioidic interface
961 if (oppAngle - angle < 0)
962 {
963 if (ExactOppMove < Exactdiff && tmp3 != tmp1)
964 {
965 RotateAngle =
966 -SectorAngle + tmp2 * SectorAngle - tmp1 * SectorAngle;
967
968 zone->Rotate(gloCoord, RotateAngle);
969 }
970 else
971 {
972 RotateAngle = -tmp1 * SectorAngle + tmp2 * SectorAngle;
973 zone->Rotate(gloCoord, RotateAngle);
974 }
975 }
976 else
977 {
978 if (Exactdiff < ExactOppMove && tmp3 == tmp1)
979 {
980 RotateAngle =
981 SectorAngle + tmp2 * SectorAngle - tmp1 * SectorAngle;
982 zone->Rotate(gloCoord, RotateAngle);
983 }
984 else
985 {
986 RotateAngle = -tmp1 * SectorAngle + tmp2 * SectorAngle;
987 zone->Rotate(gloCoord, RotateAngle);
988 }
989 }
990 }
991}
const Array< OneD, NekDouble > base(3, 0.0)

References Nektar::ErrorUtil::efatal, Nektar::SpatialDomains::eTranslate, m_interface, m_movement, and NEKERROR.

Referenced by CalcLocalMissing().

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

549{
550 // If flagged then fill trace from local coords
551 if (m_checkLocal)
552 {
553 for (auto &foundLocCoord : m_foundLocalCoords)
554 {
555 int traceId = m_trace->GetElmtToExpId(foundLocCoord.second.first);
556 Array<OneD, NekDouble> locCoord = foundLocCoord.second.second;
557
558 Array<OneD, NekDouble> edgePhys =
559 Fwd + m_trace->GetPhys_Offset(traceId);
560
561 Bwd[foundLocCoord.first] =
562 m_trace->GetExp(traceId)->StdPhysEvaluate(locCoord, edgePhys);
563 }
564 }
565}

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

688{
689 for (int i = 0; i < m_mapMissingCoordToTrace.size(); ++i)
690 {
691 if (!std::isnan(trace[i]))
692 {
693 Bwd[m_mapMissingCoordToTrace[i]] = trace[i];
694 }
695 }
696}

References m_mapMissingCoordToTrace.

◆ GetInterface()

SpatialDomains::InterfaceShPtr Nektar::MultiRegions::InterfaceTrace::GetInterface ( )
inline

Returns the interface object.

Definition at line 84 of file InterfaceMapDG.h.

85 {
86 return m_interface;
87 }

References m_interface.

◆ GetLocalCoords()

std::map< int, std::pair< int, Array< OneD, NekDouble > > > Nektar::MultiRegions::InterfaceTrace::GetLocalCoords ( )
inline

Definition at line 73 of file InterfaceMapDG.h.

74 {
75 return m_foundLocalCoords;
76 }

References m_foundLocalCoords.

◆ GetMissingCoords()

std::vector< Array< OneD, NekDouble > > Nektar::MultiRegions::InterfaceTrace::GetMissingCoords ( )
inline

Returns the missing coordinates vector.

Definition at line 67 of file InterfaceMapDG.h.

68 {
69 return m_missingCoords;
70 }

References m_missingCoords.

◆ GetRotAngle()

std::map< int, std::pair< int, NekDouble > > Nektar::MultiRegions::InterfaceTrace::GetRotAngle ( )
inline

Definition at line 78 of file InterfaceMapDG.h.

79 {
80 return m_phyRotAngle;
81 }

References m_phyRotAngle.

◆ RotLocalBwdDeriveTrace()

void Nektar::MultiRegions::InterfaceTrace::RotLocalBwdDeriveTrace ( TensorOfArray3D< NekDouble > &  Bwd,
const int &  dim 
)

Rotates the interface derivative velocity for sector rotation.

Definition at line 603 of file InterfaceMapDG.cpp.

605{
606
607 // Get rotating information from zones
608 size_t nVariables = Bwd[0].size();
609 int id = m_interface->GetId();
610 auto zones = m_movement->GetZones();
611 auto zone = std::static_pointer_cast<SpatialDomains::ZoneRotate>(zones[id]);
612 Array<OneD, NekDouble> tmpVelocity(dim);
613 auto axis = zone->GetAxis();
614 int dir = 0;
615 if (axis[0] == 0 && axis[1] == 0 && axis[2] == 1)
616 {
617 dir = 2;
618 }
619 else if (axis[0] == 1 && axis[1] == 0 && axis[2] == 0)
620 {
621 dir = 0;
622 }
623 else if (axis[0] == 0 && axis[1] == 1 && axis[2] == 0)
624 {
625 dir = 1;
626 }
627 else
628 {
629 NEKERROR(ErrorUtil::efatal, "Invalid rotation axis.");
630 }
631
632 // Rotates the derivative velocity for sector rotation
633 for (auto &phyRotAngle : m_phyRotAngle)
634 {
635 int traceId = phyRotAngle.first;
636 NekDouble angle = m_phyRotAngle[traceId].second;
637 if (angle != 0)
638 {
639 TensorOfArray3D<NekDouble> tmpBwd(dim);
640
641 for (int j = 0; j < dim; ++j)
642 {
643 tmpBwd[j] = Array<OneD, Array<OneD, NekDouble>>{nVariables};
644 for (int i = 0; i < nVariables; ++i)
645 {
646 tmpBwd[j][i] = Array<OneD, NekDouble>{1, 0.0};
647 }
648 }
649
650 for (int j = 0; j < dim; ++j)
651 {
652 for (int i = 0; i < nVariables; ++i)
653 {
654 tmpBwd[j][i][0] = Bwd[j][i][traceId];
655 }
656 }
657
658 DeriveRotate(tmpBwd, dir, -angle);
659
660 for (int j = 0; j < dim; ++j)
661 {
662 for (int i = 0; i < dim; ++i)
663 {
664 Bwd[j][i][traceId] = tmpBwd[j][i][0];
665 }
666 }
667 }
668 }
669}
void DeriveRotate(TensorOfArray3D< NekDouble > &Bwd, const int dir, const NekDouble angle)
Rotate the 3D vector field in Bwd around the axis 'dir' by 'angle'.

References DeriveRotate(), Nektar::ErrorUtil::efatal, m_interface, m_movement, m_phyRotAngle, and NEKERROR.

◆ RotLocalBwdTrace()

void Nektar::MultiRegions::InterfaceTrace::RotLocalBwdTrace ( Array< OneD, Array< OneD, NekDouble > > &  Bwd,
const int &  dim 
)

Rotates the interface velocity for sector rotation.

Definition at line 567 of file InterfaceMapDG.cpp.

569{
570
571 Array<OneD, Array<OneD, NekDouble>> vecLocs(1);
572 vecLocs[0] = Array<OneD, NekDouble>(dim);
573 for (int i = 0; i < dim; ++i)
574 {
575 vecLocs[0][i] = 1 + i;
576 }
577 int id = m_interface->GetId();
578 auto zones = m_movement->GetZones();
579 auto zone = std::static_pointer_cast<SpatialDomains::ZoneRotate>(zones[id]);
580 Array<OneD, NekDouble> tmpVelocity(dim);
581
582 // Rotates the interface velocity for sector rotation
583 // Loop over rotated search coordinates
584 for (auto &phyRotAngle : m_phyRotAngle)
585 {
586 int traceId = phyRotAngle.first;
587 NekDouble angle = m_phyRotAngle[traceId].second;
588 if (angle != 0)
589 {
590 for (int i = 0; i < dim; ++i)
591 {
592 tmpVelocity[i] = Bwd[vecLocs[0][i]][traceId];
593 }
594 zone->Rotate(tmpVelocity, -angle);
595 for (int i = 0; i < dim; ++i)
596 {
597 Bwd[vecLocs[0][i]][traceId] = tmpVelocity[i];
598 }
599 }
600 }
601}

References m_interface, m_movement, and m_phyRotAngle.

◆ SetCheckLocal()

void Nektar::MultiRegions::InterfaceTrace::SetCheckLocal ( bool  flag)
inline

Definition at line 61 of file InterfaceMapDG.h.

62 {
63 m_checkLocal = flag;
64 }

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

Referenced by CalcLocalMissing(), FillLocalBwdTrace(), and GetLocalCoords().

◆ m_interface

SpatialDomains::InterfaceShPtr Nektar::MultiRegions::InterfaceTrace::m_interface
private

Local interface object.

Definition at line 113 of file InterfaceMapDG.h.

Referenced by CalcLocalMissing(), DomainCheck(), GetInterface(), RotLocalBwdDeriveTrace(), and RotLocalBwdTrace().

◆ m_mapMissingCoordToTrace

std::vector<int> Nektar::MultiRegions::InterfaceTrace::m_mapMissingCoordToTrace
private

Vector of indices corresponding to m_missingCoord locations in trace.

Definition at line 125 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 119 of file InterfaceMapDG.h.

Referenced by CalcLocalMissing(), and GetMissingCoords().

◆ m_movement

SpatialDomains::MovementSharedPtr Nektar::MultiRegions::InterfaceTrace::m_movement
private

Movement object associated with the non-conformal interfaces.

Definition at line 115 of file InterfaceMapDG.h.

Referenced by CalcLocalMissing(), DomainCheck(), RotLocalBwdDeriveTrace(), and RotLocalBwdTrace().

◆ m_phyRotAngle

std::map<int, std::pair<int, NekDouble> > Nektar::MultiRegions::InterfaceTrace::m_phyRotAngle
private

Map of phy rotate angle for sector rotate interfaces.

Definition at line 123 of file InterfaceMapDG.h.

Referenced by CalcLocalMissing(), GetRotAngle(), RotLocalBwdDeriveTrace(), and RotLocalBwdTrace().

◆ m_trace

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

Trace expansion list.

Definition at line 111 of file InterfaceMapDG.h.

Referenced by CalcLocalMissing(), and FillLocalBwdTrace().