Nektar++
Public Member Functions | Static Public Member Functions | Static Public Attributes | Protected Member Functions | Protected Attributes | Private Member Functions | List of all members
Nektar::SolverUtils::CouplingCwipi Class Reference

#include <CouplingCwipi.h>

Inheritance diagram for Nektar::SolverUtils::CouplingCwipi:
[legend]

Public Member Functions

SOLVER_UTILS_EXPORT CouplingCwipi (MultiRegions::ExpListSharedPtr field)
 
SOLVER_UTILS_EXPORT ~CouplingCwipi () override
 
SOLVER_UTILS_EXPORT void SendCallback (Array< OneD, Array< OneD, NekDouble > > &interpField, Array< OneD, Array< OneD, NekDouble > > &distCoords)
 
- Public Member Functions inherited from Nektar::SolverUtils::Coupling
virtual SOLVER_UTILS_EXPORT ~Coupling ()
 
SOLVER_UTILS_EXPORT void Init ()
 
SOLVER_UTILS_EXPORT const std::map< std::string, std::string > GetConfig ()
 
SOLVER_UTILS_EXPORT std::vector< std::string > GetSendFieldNames ()
 
SOLVER_UTILS_EXPORT std::vector< std::string > GetRecvFieldNames ()
 
SOLVER_UTILS_EXPORT void Finalize ()
 
SOLVER_UTILS_EXPORT void Send (const int step, const NekDouble time, const Array< OneD, const Array< OneD, NekDouble > > &field, std::vector< std::string > &varNames)
 
SOLVER_UTILS_EXPORT void Receive (const int step, const NekDouble time, Array< OneD, Array< OneD, NekDouble > > &field, std::vector< std::string > &varNames)
 

Static Public Member Functions

static CouplingSharedPtr create (MultiRegions::ExpListSharedPtr field)
 Creates an instance of this class. More...
 
static SOLVER_UTILS_EXPORT void InterpCallback (const int entities_dim, const int n_local_vertex, const int n_local_element, const int n_local_polhyedra, const int n_distant_point, const double local_coordinates[], const int local_connectivity_index[], const int local_connectivity[], const int local_polyhedra_face_index[], const int local_polyhedra_cell_to_face_connectivity[], const int local_polyhedra_face_connectivity_index[], const int local_polyhedra_face_connectivity[], const double distant_points_coordinates[], const int distant_points_location[], const float distant_points_distance[], const int distant_points_barycentric_coordinates_index[], const double distant_points_barycentric_coordinates[], const int stride, const cwipi_solver_type_t solver_type, const void *local_field, void *distant_field)
 

Static Public Attributes

static std::string className
 

Protected Member Functions

SOLVER_UTILS_EXPORT void v_Init () override
 
SOLVER_UTILS_EXPORT void v_Send (const int step, const NekDouble time, const Array< OneD, const Array< OneD, NekDouble > > &field, std::vector< std::string > &varNames) override
 
SOLVER_UTILS_EXPORT void v_Receive (const int step, const NekDouble time, Array< OneD, Array< OneD, NekDouble > > &field, std::vector< std::string > &varNames) override
 
SOLVER_UTILS_EXPORT void v_Finalize () override
 
SOLVER_UTILS_EXPORT NekDouble GetSendField (const int i, const int j) const
 
- Protected Member Functions inherited from Nektar::SolverUtils::Coupling
SOLVER_UTILS_EXPORT Coupling (MultiRegions::ExpListSharedPtr field)
 
virtual SOLVER_UTILS_EXPORT void v_Init ()
 
virtual SOLVER_UTILS_EXPORT void v_Send (const int step, const NekDouble time, const Array< OneD, const Array< OneD, NekDouble > > &field, std::vector< std::string > &varNames)=0
 
virtual SOLVER_UTILS_EXPORT void v_Receive (const int step, const NekDouble time, Array< OneD, Array< OneD, NekDouble > > &field, std::vector< std::string > &varNames)=0
 
virtual SOLVER_UTILS_EXPORT void v_Finalize ()
 
SOLVER_UTILS_EXPORT std::vector< int > GenerateVariableMapping (std::vector< std::string > &vars, std::vector< std::string > &transVars)
 

Protected Attributes

NekDouble m_filtWidth
 
Array< OneD, Array< OneD, NekDouble > > m_sendField
 
MultiRegions::ExpListSharedPtr m_recvField
 
Array< OneD, Array< OneD, NekDouble > > m_oldFields
 
Array< OneD, Array< OneD, NekDouble > > m_newFields
 
int m_sendHandle
 
int m_recvHandle
 
int m_sendTag
 
int m_recvTag
 
int m_lastSend
 
int m_lastReceive
 
int m_spacedim
 
int m_nPoints
 
double * m_points
 
double * m_coords
 
int * m_connecIdx
 
int * m_connec
 
double * m_rValsInterl
 
double * m_sValsInterl
 
std::map< int, int > m_vertMap
 
std::shared_ptr< FieldUtils::Interpolator< std::vector< MultiRegions::ExpListSharedPtr > > > m_sendInterpolator
 
std::shared_ptr< FieldUtils::Interpolator< std::vector< MultiRegions::ExpListSharedPtr > > > m_extrapInterpolator
 
- Protected Attributes inherited from Nektar::SolverUtils::Coupling
std::string m_couplingName
 
CouplingConfigMap m_config
 
MultiRegions::ExpListSharedPtr m_evalField
 
int m_nSendVars
 
std::vector< std::string > m_sendFieldNames
 
int m_sendSteps
 
int m_nRecvVars
 
std::vector< std::string > m_recvFieldNames
 
int m_recvSteps
 

Private Member Functions

void ReadConfig (LibUtilities::SessionReaderSharedPtr session)
 
void SetupReceive ()
 
void SetupSend ()
 
void SendComplete ()
 
void ReceiveStart ()
 
void ReceiveCwipi (const int step, const NekDouble time, Array< OneD, Array< OneD, NekDouble > > &field)
 
void EvaluateFields (Array< OneD, Array< OneD, NekDouble > > interpField, Array< OneD, Array< OneD, NekDouble > > distCoords)
 
void SetupSendInterpolation ()
 
void AnnounceMesh ()
 
void DumpRawFields (const NekDouble time, Array< OneD, Array< OneD, NekDouble > > rVals)
 
void ExtrapolateFields (Array< OneD, Array< OneD, NekDouble > > &rVals, Array< OneD, int > &notLoc)
 
template<typename T >
void AddElementsToMesh (T geom, int &coordsPos, int &connecPos, int &conidxPos)
 

Additional Inherited Members

- Public Types inherited from Nektar::SolverUtils::Coupling
typedef std::map< std::string, std::string > CouplingConfigMap
 

Detailed Description

Definition at line 49 of file CouplingCwipi.h.

Constructor & Destructor Documentation

◆ CouplingCwipi()

Nektar::SolverUtils::CouplingCwipi::CouplingCwipi ( MultiRegions::ExpListSharedPtr  field)

Definition at line 106 of file CouplingCwipi.cpp.

107 : Coupling(field), m_sendHandle(-1), m_recvHandle(-1), m_lastSend(-1E6),
108 m_lastReceive(-1E6), m_points(nullptr), m_coords(nullptr),
109 m_connecIdx(nullptr), m_connec(nullptr), m_rValsInterl(nullptr),
110 m_sValsInterl(nullptr)
111{
112 // defaults
113 m_config["GEOMTOL"] = "0.1";
114 m_config["LOCALNAME"] = "nektar";
115 m_config["REMOTENAME"] = "precise";
116 m_config["OVERSAMPLE"] = "0";
117 m_config["FILTERWIDTH"] = "-1";
118 m_config["DUMPRAW"] = "0";
119 m_config["SENDMETHOD"] = "EVALUATE";
120 m_config["NOTLOCMETHOD"] = "KEEP";
121}
CouplingConfigMap m_config
Definition: Coupling.h:109
SOLVER_UTILS_EXPORT Coupling(MultiRegions::ExpListSharedPtr field)
Definition: Coupling.cpp:48

References Nektar::SolverUtils::Coupling::m_config.

◆ ~CouplingCwipi()

Nektar::SolverUtils::CouplingCwipi::~CouplingCwipi ( )
override

Definition at line 213 of file CouplingCwipi.cpp.

214{
215 free(m_coords);
216 free(m_points);
217 free(m_connec);
218 free(m_connecIdx);
219 free(m_rValsInterl);
220 free(m_sValsInterl);
221}

References m_connec, m_connecIdx, m_coords, m_points, m_rValsInterl, and m_sValsInterl.

Member Function Documentation

◆ AddElementsToMesh()

template<typename T >
void Nektar::SolverUtils::CouplingCwipi::AddElementsToMesh ( geom,
int &  coordsPos,
int &  connecPos,
int &  conidxPos 
)
private

Definition at line 495 of file CouplingCwipi.cpp.

497{
498 // helper variables
499 Array<OneD, NekDouble> x(3);
501 int vertID;
502
503 int kNverts = T::mapped_type::element_type::kNverts;
504
505 // iterate over all elements
506 for (auto it = geom.begin(); it != geom.end(); it++)
507 {
508 // iterate over the elements vertices
509 for (int j = 0; j < kNverts; ++j)
510 {
511 vert = it->second->GetVertex(j);
512 vertID = vert->GetVid();
513
514 // check if we already stored the vertex
515 if (m_vertMap.count(vertID) == 0)
516 {
517 // store the vertex
518 vert->GetCoords(x[0], x[1], x[2]);
519 m_coords[3 * coordsPos + 0] = double(x[0]);
520 m_coords[3 * coordsPos + 1] = double(x[1]);
521 m_coords[3 * coordsPos + 2] = double(x[2]);
522
523 // store the vertex position in the m_coords array
524 m_vertMap[vertID] = coordsPos;
525 coordsPos++;
526 }
527
528 m_connec[connecPos] = m_vertMap[vertID] + 1;
529 connecPos++;
530 }
531
532 m_connecIdx[conidxPos + 1] = m_connecIdx[conidxPos] + kNverts;
533 conidxPos++;
534 }
535}
std::shared_ptr< PointGeom > PointGeomSharedPtr
Definition: Geometry.h:57

References m_connec, m_connecIdx, m_coords, and m_vertMap.

Referenced by AnnounceMesh().

◆ AnnounceMesh()

void Nektar::SolverUtils::CouplingCwipi::AnnounceMesh ( )
private

Definition at line 409 of file CouplingCwipi.cpp.

410{
412
413 // get Elements
421 if (m_spacedim == 1)
422 {
423 seggeom = graph->GetAllSegGeoms();
424 }
425 else if (m_spacedim == 2)
426 {
427 trigeom = graph->GetAllTriGeoms();
428 quadgeom = graph->GetAllQuadGeoms();
429 }
430 else if (m_spacedim == 3)
431 {
432 tetgeom = graph->GetAllTetGeoms();
433 pyrgeom = graph->GetAllPyrGeoms();
434 prismgeom = graph->GetAllPrismGeoms();
435 hexgeom = graph->GetAllHexGeoms();
436 };
437
438 int nVerts = graph->GetNvertices();
439 int nElts = seggeom.size() + trigeom.size() + quadgeom.size() +
440 tetgeom.size() + pyrgeom.size() + prismgeom.size() +
441 hexgeom.size();
442
443 // allocate CWIPI arrays
444 m_coords = (double *)malloc(sizeof(double) * 3 * nVerts);
445 ASSERTL1(m_coords != nullptr, "malloc failed for m_coords");
446 int tmp = 2 * seggeom.size() + 3 * trigeom.size() + 4 * quadgeom.size() +
447 4 * tetgeom.size() + 5 * pyrgeom.size() + 6 * prismgeom.size() +
448 8 * hexgeom.size();
449 m_connec = (int *)malloc(sizeof(int) * tmp);
450 ASSERTL1(m_connec != nullptr, "malloc failed for m_connec");
451 m_connecIdx = (int *)malloc(sizeof(int) * (nElts + 1));
452 ASSERTL1(m_connecIdx != nullptr, "malloc failed for m_connecIdx");
453
454 m_connecIdx[0] = 0;
455 int coordsPos = 0;
456 int connecPos = 0;
457 int conidxPos = 0;
458
459 AddElementsToMesh(seggeom, coordsPos, connecPos, conidxPos);
460 AddElementsToMesh(trigeom, coordsPos, connecPos, conidxPos);
461 AddElementsToMesh(quadgeom, coordsPos, connecPos, conidxPos);
462 AddElementsToMesh(tetgeom, coordsPos, connecPos, conidxPos);
463 AddElementsToMesh(pyrgeom, coordsPos, connecPos, conidxPos);
464 AddElementsToMesh(prismgeom, coordsPos, connecPos, conidxPos);
465 AddElementsToMesh(hexgeom, coordsPos, connecPos, conidxPos);
466
467 // output the mesh in tecplot format. If this works, CWIPI will be able
468 // to process it, too
469 /*
470 cout << "VARIABLES = \"X\", \"Y\", \"Z\", \"U\"" << endl;
471 cout << "ZONE NODES=" << nVerts << ",ELEMENTS=" << nElts;
472 cout << ",DATAPACKING=POINT, ZONETYPE=FETETRAHEDRON" << endl;
473 for (int i = 0; i < nVerts; ++i)
474 {
475 cout << m_coords[3*i + 0] << " " << m_coords[3*i + 1] << " ";
476 cout << m_coords[3*i + 2] << " " << 1.0 << endl;
477 }
478 for (int i = 0; i < nElts; ++i)
479 {
480 cout << m_connec[i*4 + 0] << " " << m_connec[i*4 + 1] << " ";
481 cout << m_connec[i*4 + 2] << " " << m_connec[i*4 + 3] << endl;
482 }
483 */
484
485 cwipi_define_mesh(m_couplingName.c_str(), nVerts, nElts, m_coords,
487}
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:242
void AddElementsToMesh(T geom, int &coordsPos, int &connecPos, int &conidxPos)
MultiRegions::ExpListSharedPtr m_evalField
Definition: Coupling.h:111
std::map< int, TriGeomSharedPtr > TriGeomMap
Definition: TriGeom.h:57
std::map< int, PyrGeomSharedPtr > PyrGeomMap
Definition: PyrGeom.h:76
std::map< int, QuadGeomSharedPtr > QuadGeomMap
Definition: QuadGeom.h:52
std::map< int, SegGeomSharedPtr > SegGeomMap
Definition: SegGeom.h:49
std::map< int, TetGeomSharedPtr > TetGeomMap
Definition: TetGeom.h:84
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:174
std::map< int, PrismGeomSharedPtr > PrismGeomMap
Definition: PrismGeom.h:83
std::map< int, HexGeomSharedPtr > HexGeomMap
Definition: HexGeom.h:85

References AddElementsToMesh(), ASSERTL1, m_connec, m_connecIdx, m_coords, Nektar::SolverUtils::Coupling::m_couplingName, Nektar::SolverUtils::Coupling::m_evalField, and m_spacedim.

Referenced by v_Init().

◆ create()

static CouplingSharedPtr Nektar::SolverUtils::CouplingCwipi::create ( MultiRegions::ExpListSharedPtr  field)
inlinestatic

Creates an instance of this class.

Definition at line 56 of file CouplingCwipi.h.

57 {
60 p->Init();
61 return p;
62 }
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
SOLVER_UTILS_EXPORT typedef std::shared_ptr< Coupling > CouplingSharedPtr
Definition: Coupling.h:45

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::SolverUtils::CouplingSharedPtr, and CellMLToNektar.cellml_metadata::p.

◆ DumpRawFields()

void Nektar::SolverUtils::CouplingCwipi::DumpRawFields ( const NekDouble  time,
Array< OneD, Array< OneD, NekDouble > >  rVals 
)
private

Definition at line 1019 of file CouplingCwipi.cpp.

1021{
1022 LibUtilities::Timer timer1;
1023 timer1.Start();
1024
1025#if (defined _WIN32 && _MSC_VER < 1900)
1026 // We need this to make sure boost::format has always
1027 // two digits in the exponents of Scientific notation.
1028 unsigned int old_exponent_format;
1029 old_exponent_format = _set_output_format(_TWO_DIGIT_EXPONENT);
1030 filename = boost::str(boost::format(filename) % m_time);
1031 _set_output_format(old_exponent_format);
1032#else
1033 std::string filename =
1034 boost::str(boost::format(m_config["DUMPRAW"]) % time);
1035#endif
1036
1037 Array<OneD, Array<OneD, NekDouble>> tmp(m_nRecvVars + m_spacedim);
1038 for (int i = 0; i < 3; ++i)
1039 {
1040 tmp[i] = Array<OneD, NekDouble>(m_recvField->GetTotPoints(), 0.0);
1041 }
1042 m_recvField->GetCoords(tmp[0], tmp[1], tmp[2]);
1043
1044 for (int i = 0; i < m_nRecvVars; ++i)
1045 {
1046 tmp[m_spacedim + i] = rVals[i];
1047 }
1048
1049 LibUtilities::CsvIO csvIO(m_evalField->GetSession()->GetComm());
1053 csvIO.Write(filename, rvPts);
1054
1055 timer1.Stop();
1056 if (m_evalField->GetComm()->GetRank() == 0 &&
1057 m_evalField->GetSession()->DefinesCmdLineArgument("verbose"))
1058 {
1059 cout << "DumpRawFields total time: " << timer1.TimePerTest(1) << endl;
1060 }
1061}
MultiRegions::ExpListSharedPtr m_recvField
Definition: CouplingCwipi.h:94
std::vector< std::string > m_recvFieldNames
Definition: Coupling.h:118
std::shared_ptr< PtsField > PtsFieldSharedPtr
Definition: PtsField.h:184

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), CellMLToNektar.pycml::format, Nektar::SolverUtils::Coupling::m_config, Nektar::SolverUtils::Coupling::m_evalField, Nektar::SolverUtils::Coupling::m_nRecvVars, m_recvField, Nektar::SolverUtils::Coupling::m_recvFieldNames, m_spacedim, Nektar::LibUtilities::Timer::Start(), Nektar::LibUtilities::Timer::Stop(), Nektar::LibUtilities::Timer::TimePerTest(), and Nektar::LibUtilities::CsvIO::Write().

Referenced by ReceiveCwipi().

◆ EvaluateFields()

void Nektar::SolverUtils::CouplingCwipi::EvaluateFields ( Array< OneD, Array< OneD, NekDouble > >  interpField,
Array< OneD, Array< OneD, NekDouble > >  distCoords 
)
private

Definition at line 320 of file CouplingCwipi.cpp.

323{
324 int nOutPts = distCoords.size();
325
326 Array<OneD, NekDouble> Lcoords(m_spacedim, 0.0);
327 for (int i = 0; i < nOutPts; ++i)
328 {
329 // Obtain Element and LocalCoordinate to interpolate
330 int elmtid = -1;
332 while (elmtid < 0 && tol <= 1E3 * NekConstants::kNekZeroTol)
333 {
334 elmtid = m_evalField->GetExpIndex(distCoords[i], Lcoords, tol);
335 tol *= 2;
336 }
337 if (tol > 2 * NekConstants::kNekZeroTol)
338 {
339 for (int j = 0; j < m_spacedim; ++j)
340 {
341 if (Lcoords[j] < -1 - 0.75 * NekConstants::kNekZeroTol)
342 {
343 Lcoords[j] = -1;
344 }
345 if (Lcoords[j] > 1 + 0.75 * NekConstants::kNekZeroTol)
346 {
347 Lcoords[j] = 1;
348 }
349 }
350 }
351
352 ASSERTL0(elmtid >= 0,
353 "no element found for (" +
354 boost::lexical_cast<string>(distCoords[i][0]) + ", " +
355 boost::lexical_cast<string>(distCoords[i][1]) + ", " +
356 boost::lexical_cast<string>(distCoords[i][2]) + ")");
357
358 int offset = m_evalField->GetPhys_Offset(elmtid);
359
360 for (int f = 0; f < m_nSendVars; ++f)
361 {
362 NekDouble value = m_evalField->GetExp(elmtid)->StdPhysEvaluate(
363 Lcoords, m_sendField[f] + offset);
364
365 ASSERTL0(!(boost::math::isnan)(value), "new value is not a number");
366 interpField[f][i] = value;
367 }
368 }
369}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
Array< OneD, Array< OneD, NekDouble > > m_sendField
Definition: CouplingCwipi.h:92
static const NekDouble kNekZeroTol
double NekDouble

References ASSERTL0, Nektar::NekConstants::kNekZeroTol, Nektar::SolverUtils::Coupling::m_evalField, Nektar::SolverUtils::Coupling::m_nSendVars, m_sendField, and m_spacedim.

Referenced by SendCallback().

◆ ExtrapolateFields()

void Nektar::SolverUtils::CouplingCwipi::ExtrapolateFields ( Array< OneD, Array< OneD, NekDouble > > &  rVals,
Array< OneD, int > &  notLoc 
)
private

Definition at line 893 of file CouplingCwipi.cpp.

895{
896 LibUtilities::Timer timer1, timer2;
897 timer1.Start();
898
899 int totvars = 3 + m_nRecvVars;
900 int nNotLoc = notLoc.size();
901 int nranks = m_evalField->GetSession()->GetComm()->GetSize();
902
903 Array<OneD, int> thisNLoc(1, m_nPoints - nNotLoc);
904 Array<OneD, int> sizeMap(nranks);
905 m_evalField->GetSession()->GetComm()->AllGather(thisNLoc, sizeMap);
906
907 Array<OneD, int> offsetMap(nranks);
908 offsetMap[0] = 0;
909 for (int i = 1; i < nranks; ++i)
910 {
911 offsetMap[i] = offsetMap[i - 1] + sizeMap[i - 1];
912 }
913 int totNLoc = offsetMap[nranks - 1] + sizeMap[nranks - 1];
914
915 Array<OneD, Array<OneD, NekDouble>> allVals(totvars);
916 for (int i = 0; i < 3; ++i)
917 {
918 allVals[i] = Array<OneD, NekDouble>(m_nPoints);
919 }
920 m_recvField->GetCoords(allVals[0], allVals[1], allVals[2]);
921
922 if (m_spacedim < 3)
923 {
924 Vmath::Zero(m_nPoints, allVals[2], 1);
925 }
926 if (m_spacedim < 2)
927 {
928 Vmath::Zero(m_nPoints, allVals[1], 1);
929 }
930
931 for (int i = 0; i < m_nRecvVars; ++i)
932 {
933 allVals[3 + i] = rVals[i];
934 }
935
936 Array<OneD, Array<OneD, NekDouble>> locatedVals(totvars);
937 for (int i = 0; i < totvars; ++i)
938 {
939 locatedVals[i] = Array<OneD, NekDouble>(totNLoc, -42.0);
940 }
941
942 // only copy points from allVals to locatedVals that were located
943 int offset = offsetMap[m_evalField->GetSession()->GetComm()->GetRank()];
944 for (int i = 0, intPos = 0, locPos = 0; i < m_nPoints; ++i)
945 {
946 if (locPos < nNotLoc && notLoc[locPos] == i)
947 {
948 // do nothing
949 locPos++;
950 }
951 else
952 {
953 for (int k = 0; k < totvars; ++k)
954 {
955 locatedVals[k][offset + intPos] = allVals[k][i];
956 }
957 intPos++;
958 }
959 }
960
961 // send all located points to all ranks. This is probably horribly expensive
962 timer2.Start();
963 for (int i = 0; i < totvars; ++i)
964 {
965 m_evalField->GetSession()->GetComm()->AllGatherv(locatedVals[i],
966 sizeMap, offsetMap);
967 }
968 timer2.Stop();
969
970 if (nNotLoc > 0)
971 {
974 3, locatedVals);
975
976 Array<OneD, Array<OneD, NekDouble>> notLocVals(totvars);
977 for (int j = 0; j < totvars; ++j)
978 {
979 notLocVals[j] = Array<OneD, NekDouble>(nNotLoc);
980 for (int i = 0; i < nNotLoc; ++i)
981 {
982 notLocVals[j][i] = allVals[j][notLoc[i]];
983 }
984 }
987 3, notLocVals);
988
989 // perform a nearest neighbour interpolation from locatedVals to the not
990 // located rVals
992 {
993 m_extrapInterpolator = MemoryManager<FieldUtils::Interpolator<
994 std::vector<MultiRegions::ExpListSharedPtr>>>::
995 AllocateSharedPtr(LibUtilities::eNearestNeighbour);
996 m_extrapInterpolator->CalcWeights(locatedPts, notlocPts);
997 m_extrapInterpolator->PrintStatistics();
998 }
999 m_extrapInterpolator->Interpolate(locatedPts, notlocPts);
1000
1001 for (int j = 3; j < totvars; ++j)
1002 {
1003 for (int i = 0; i < nNotLoc; ++i)
1004 {
1005 allVals[j][notLoc[i]] = notLocVals[j][i];
1006 }
1007 }
1008 }
1009
1010 timer1.Stop();
1011 if (m_evalField->GetComm()->GetRank() == 0 &&
1012 m_evalField->GetSession()->DefinesCmdLineArgument("verbose"))
1013 {
1014 cout << "ExtrapolateFields total time: " << timer1.TimePerTest(1);
1015 cout << " (AllGathervI: " << timer2.TimePerTest(1) << ")" << endl;
1016 }
1017}
std::shared_ptr< FieldUtils::Interpolator< std::vector< MultiRegions::ExpListSharedPtr > > > m_extrapInterpolator
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.hpp:273

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::LibUtilities::eNearestNeighbour, Nektar::SolverUtils::Coupling::m_evalField, m_extrapInterpolator, m_nPoints, Nektar::SolverUtils::Coupling::m_nRecvVars, m_recvField, m_spacedim, Nektar::LibUtilities::Timer::Start(), Nektar::LibUtilities::Timer::Stop(), Nektar::LibUtilities::Timer::TimePerTest(), and Vmath::Zero().

Referenced by ReceiveCwipi().

◆ GetSendField()

SOLVER_UTILS_EXPORT NekDouble Nektar::SolverUtils::CouplingCwipi::GetSendField ( const int  i,
const int  j 
) const
inlineprotected

Definition at line 141 of file CouplingCwipi.h.

142 {
143 return m_sendField[i][j];
144 }

References m_sendField.

◆ InterpCallback()

void Nektar::SolverUtils::CouplingCwipi::InterpCallback ( const int  entities_dim,
const int  n_local_vertex,
const int  n_local_element,
const int  n_local_polhyedra,
const int  n_distant_point,
const double  local_coordinates[],
const int  local_connectivity_index[],
const int  local_connectivity[],
const int  local_polyhedra_face_index[],
const int  local_polyhedra_cell_to_face_connectivity[],
const int  local_polyhedra_face_connectivity_index[],
const int  local_polyhedra_face_connectivity[],
const double  distant_points_coordinates[],
const int  distant_points_location[],
const float  distant_points_distance[],
const int  distant_points_barycentric_coordinates_index[],
const double  distant_points_barycentric_coordinates[],
const int  stride,
const cwipi_solver_type_t  solver_type,
const void local_field,
void distant_field 
)
static

Definition at line 59 of file CouplingCwipi.cpp.

77{
78 Array<OneD, Array<OneD, NekDouble>> interpField(stride);
79
80 Array<OneD, Array<OneD, NekDouble>> distCoords(n_distant_point);
81 for (int i = 0; i < n_distant_point; ++i)
82 {
83 distCoords[i] = Array<OneD, NekDouble>(3);
84 for (int j = 0; j < 3; ++j)
85 {
86 distCoords[i][j] = distant_points_coordinates[3 * i + j];
87 }
88 }
89
90 std::stringstream sst;
91 sst << entities_dim << "," << n_local_vertex << "," << stride;
92 SendCallbackMap[sst.str()](interpField, distCoords);
93
94 ASSERTL0(interpField.size() == stride, "size mismatch");
95 ASSERTL0(interpField[0].size() == n_distant_point, "size mismatch");
96
97 for (int i = 0; i < n_distant_point; i++)
98 {
99 for (int j = 0; j < stride; ++j)
100 {
101 ((double *)distant_field)[i * stride + j] = interpField[j][i];
102 }
103 }
104}
static std::map< std::string, SendCallbackType > SendCallbackMap

References ASSERTL0, and Nektar::SolverUtils::SendCallbackMap.

Referenced by SetupSend().

◆ ReadConfig()

void Nektar::SolverUtils::CouplingCwipi::ReadConfig ( LibUtilities::SessionReaderSharedPtr  session)
private

Definition at line 223 of file CouplingCwipi.cpp.

224{
225 ASSERTL0(session->DefinesElement("Nektar/Coupling"),
226 "No Coupling config found");
227
228 m_config["LOCALNAME"] = session->GetCmdLineArgument<std::string>("cwipi");
229
230 TiXmlElement *vCoupling = session->GetElement("Nektar/Coupling");
231 ASSERTL0(vCoupling, "Invalid Coupling config");
232
233 m_filtWidth = boost::lexical_cast<NekDouble>(m_config["FILTERWIDTH"]);
234}

References ASSERTL0, Nektar::SolverUtils::Coupling::m_config, and m_filtWidth.

Referenced by v_Init().

◆ ReceiveCwipi()

void Nektar::SolverUtils::CouplingCwipi::ReceiveCwipi ( const int  step,
const NekDouble  time,
Array< OneD, Array< OneD, NekDouble > > &  field 
)
private

Definition at line 716 of file CouplingCwipi.cpp.

718{
719 ASSERTL1(m_nRecvVars == field.size(), "field size mismatch");
720
721 if (m_nRecvVars < 1 || m_recvSteps < 1)
722 {
723 return;
724 }
725
726 if (step >= m_lastReceive + m_recvSteps)
727 {
728 m_lastReceive = step;
729
730 if (m_evalField->GetComm()->GetRank() == 0 &&
731 m_evalField->GetSession()->DefinesCmdLineArgument("verbose"))
732 {
733 cout << "waiting for receive at i = " << step << ", t = " << time
734 << endl;
735 }
736
737 LibUtilities::Timer timer1, timer2, timer3;
738 timer1.Start();
739
740 Array<OneD, NekDouble> tmpC(m_recvField->GetNcoeffs());
741 Array<OneD, Array<OneD, NekDouble>> rVals(m_nRecvVars);
742 for (int i = 0; i < m_nRecvVars; ++i)
743 {
744 rVals[i] = Array<OneD, NekDouble>(m_recvField->GetTotPoints());
745 }
746
747 timer2.Start();
748 cwipi_wait_irecv(m_couplingName.c_str(), m_recvHandle);
749 timer2.Stop();
750
751 // set to -1 so we know we can start receiving again
752 m_recvHandle = -1;
753
754 int nNotLoc = cwipi_get_n_not_located_points(m_couplingName.c_str());
755 Array<OneD, int> notLoc;
756 if (nNotLoc != 0)
757 {
758 cout << "WARNING: relocating " << nNotLoc << " of " << m_nPoints
759 << " points" << endl;
760
761 const int *tmp =
762 cwipi_get_not_located_points(m_couplingName.c_str());
763 notLoc = Array<OneD, int>(nNotLoc);
764 for (int i = 0; i < nNotLoc; ++i)
765 {
766 notLoc[i] = tmp[i] - 1;
767 }
768
769 if (boost::to_upper_copy(m_config["NOTLOCMETHOD"]) == "KEEP")
770 {
771 // interpolate from m_evalField to m_recvField
772 for (int i = 0; i < m_nRecvVars; ++i)
773 {
774 m_evalField->FwdTrans(field[i], tmpC);
775 m_recvField->BwdTrans(tmpC, rVals[i]);
776 }
777 }
778 }
779
780 for (int i = 0, locPos = 0, intPos = 0; i < m_nPoints; ++i)
781 {
782 if (locPos < nNotLoc && notLoc[locPos] == i)
783 {
784 // keep the original value of field[j][i]
785 locPos++;
786 }
787 else
788 {
789 for (int j = 0; j < m_nRecvVars; ++j)
790 {
791 rVals[j][i] = m_rValsInterl[intPos * m_nRecvVars + j];
792 }
793 intPos++;
794 }
795 }
796
797 if (boost::to_upper_copy(m_config["NOTLOCMETHOD"]) == "EXTRAPOLATE")
798 {
799 int doExtrapolate = 0;
800 if (nNotLoc != 0)
801 {
802 doExtrapolate = 1;
803 }
804 m_evalField->GetSession()->GetComm()->AllReduce(
805 doExtrapolate, LibUtilities::ReduceMax);
806 if (doExtrapolate > 0)
807 {
808 ExtrapolateFields(rVals, notLoc);
809 }
810 }
811
812 if (m_config["DUMPRAW"] != "0")
813 {
814 DumpRawFields(time, rVals);
815 }
816
817 if (m_filtWidth > 0)
818 {
819 for (int i = 0; i < m_nRecvVars; ++i)
820 {
821 timer3.Start();
822
823 Array<OneD, NekDouble> forcing(m_nPoints);
824
825 Array<OneD, Array<OneD, NekDouble>> Velocity(m_spacedim);
826 for (int j = 0; j < m_spacedim; ++j)
827 {
828 Velocity[j] = Array<OneD, NekDouble>(m_nPoints, 0.0);
829 }
830
831 Vmath::Smul(m_nPoints, -m_filtWidth, rVals[i], 1, forcing, 1);
832
834 StdRegions::VarCoeffMap varcoeffs;
835 MultiRegions::VarFactorsMap varfactors =
837
839
840 // Set advection velocities
841 StdRegions::VarCoeffType varcoefftypes[] = {
844 for (int j = 0; j < m_spacedim; j++)
845 {
846 varcoeffs[varcoefftypes[j]] = Velocity[j];
847 }
848
849 // Note we are using the
850 // LinearAdvectionDiffusionReaction solver here
851 // instead of HelmSolve since m_filtWidth is negative and
852 // so matrices are not positive definite. Ideally
853 // should allow for negative m_filtWidth coefficient in
854 // HelmSolve
855 // m_recvField->LinearAdvectionDiffusionReactionSolve(
856 // Velocity, forcing, tmpC, -m_filtWidth);
857 m_recvField->LinearAdvectionDiffusionReactionSolve(
858 forcing, tmpC, factors, varcoeffs, varfactors);
859
860 m_evalField->BwdTrans(tmpC, field[i]);
861 timer3.Stop();
862
863 if (m_evalField->GetComm()->GetRank() == 0 &&
864 m_evalField->GetSession()->DefinesCmdLineArgument(
865 "verbose"))
866 {
867 cout << "Smoother time (" << m_recvFieldNames[i]
868 << "): " << timer3.TimePerTest(1) << endl;
869 }
870 }
871 }
872 else
873 {
874 for (int i = 0; i < m_nRecvVars; ++i)
875 {
876 m_recvField->FwdTrans(rVals[i], tmpC);
877 m_evalField->BwdTrans(tmpC, field[i]);
878 }
879 }
880 timer1.Stop();
881
882 if (m_evalField->GetComm()->GetRank() == 0 &&
883 m_evalField->GetSession()->DefinesCmdLineArgument("verbose"))
884 {
885 cout << "Receive total time: " << timer1.TimePerTest(1) << ", ";
886 cout << "Receive waiting time: " << timer2.TimePerTest(1) << endl;
887 }
888
889 ReceiveStart();
890 }
891}
void ExtrapolateFields(Array< OneD, Array< OneD, NekDouble > > &rVals, Array< OneD, int > &notLoc)
void DumpRawFields(const NekDouble time, Array< OneD, Array< OneD, NekDouble > > rVals)
static VarFactorsMap NullVarFactorsMap
std::map< StdRegions::ConstFactorType, Array< OneD, NekDouble > > VarFactorsMap
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:402
std::map< StdRegions::VarCoeffType, VarCoeffEntry > VarCoeffMap
Definition: StdRegions.hpp:346
StdRegions::ConstFactorMap factors
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition: Vmath.hpp:100

References ASSERTL1, DumpRawFields(), Nektar::StdRegions::eFactorLambda, Nektar::StdRegions::eVarCoeffVelX, Nektar::StdRegions::eVarCoeffVelY, Nektar::StdRegions::eVarCoeffVelZ, ExtrapolateFields(), Nektar::VarcoeffHashingTest::factors, Nektar::SolverUtils::Coupling::m_config, Nektar::SolverUtils::Coupling::m_couplingName, Nektar::SolverUtils::Coupling::m_evalField, m_filtWidth, m_lastReceive, m_nPoints, Nektar::SolverUtils::Coupling::m_nRecvVars, m_recvField, Nektar::SolverUtils::Coupling::m_recvFieldNames, m_recvHandle, Nektar::SolverUtils::Coupling::m_recvSteps, m_rValsInterl, m_spacedim, Nektar::MultiRegions::NullVarFactorsMap, ReceiveStart(), Nektar::LibUtilities::ReduceMax, Vmath::Smul(), Nektar::LibUtilities::Timer::Start(), Nektar::LibUtilities::Timer::Stop(), and Nektar::LibUtilities::Timer::TimePerTest().

Referenced by v_Receive().

◆ ReceiveStart()

void Nektar::SolverUtils::CouplingCwipi::ReceiveStart ( )
private

Definition at line 641 of file CouplingCwipi.cpp.

642{
643 if (m_recvHandle >= 0)
644 {
645 return;
646 }
647
648 LibUtilities::Timer timer1;
649 timer1.Start();
650 // workaround a bug in cwipi: receiving_field_name should be const char* but
651 // is char*
652 char recFN[10];
653 strcpy(recFN, "dummyName");
654
655 cwipi_irecv(m_couplingName.c_str(), "ex1", m_recvTag, m_nRecvVars, 0, 0.0,
656 recFN, m_rValsInterl, &m_recvHandle);
657 timer1.Stop();
658
659 if (m_evalField->GetComm()->GetRank() == 0 &&
660 m_evalField->GetSession()->DefinesCmdLineArgument("verbose"))
661 {
662 cout << "Receive start time: " << timer1.TimePerTest(1) << endl;
663 }
664}

References Nektar::SolverUtils::Coupling::m_couplingName, Nektar::SolverUtils::Coupling::m_evalField, Nektar::SolverUtils::Coupling::m_nRecvVars, m_recvHandle, m_recvTag, m_rValsInterl, Nektar::LibUtilities::Timer::Start(), Nektar::LibUtilities::Timer::Stop(), and Nektar::LibUtilities::Timer::TimePerTest().

Referenced by ReceiveCwipi(), and v_Init().

◆ SendCallback()

void Nektar::SolverUtils::CouplingCwipi::SendCallback ( Array< OneD, Array< OneD, NekDouble > > &  interpField,
Array< OneD, Array< OneD, NekDouble > > &  distCoords 
)

Definition at line 537 of file CouplingCwipi.cpp.

540{
541 ASSERTL0(interpField.size() == m_nSendVars, "size mismatch");
542
543 for (int i = 0; i < m_nSendVars; ++i)
544 {
545 interpField[i] = Array<OneD, NekDouble>(distCoords.size());
546 }
547
548 if (boost::to_upper_copy(m_config["SENDMETHOD"]) == "NEARESTNEIGHBOUR" ||
549 boost::to_upper_copy(m_config["SENDMETHOD"]) == "SHEPARD")
550 {
552 {
554 }
555
558 0, m_sendField);
561 0, interpField);
562 m_sendInterpolator->Interpolate(ptsIn, ptsOut);
563 }
564 else
565 {
566 EvaluateFields(interpField, distCoords);
567 }
568}
std::shared_ptr< FieldUtils::Interpolator< std::vector< MultiRegions::ExpListSharedPtr > > > m_sendInterpolator
void EvaluateFields(Array< OneD, Array< OneD, NekDouble > > interpField, Array< OneD, Array< OneD, NekDouble > > distCoords)

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), ASSERTL0, EvaluateFields(), Nektar::SolverUtils::Coupling::m_config, Nektar::SolverUtils::Coupling::m_nSendVars, m_sendField, m_sendInterpolator, and SetupSendInterpolation().

Referenced by SetupSend().

◆ SendComplete()

void Nektar::SolverUtils::CouplingCwipi::SendComplete ( )
private

Definition at line 619 of file CouplingCwipi.cpp.

620{
621 if (m_sendHandle < 0)
622 {
623 return;
624 }
625
626 LibUtilities::Timer timer1;
627 timer1.Start();
628 cwipi_wait_issend(m_couplingName.c_str(), m_sendHandle);
629 timer1.Stop();
630
631 // set to -1 so we dont try finishing a send before a new one was started
632 m_sendHandle = -1;
633
634 if (m_evalField->GetComm()->GetRank() == 0 &&
635 m_evalField->GetSession()->DefinesCmdLineArgument("verbose"))
636 {
637 cout << "Send waiting time: " << timer1.TimePerTest(1) << endl;
638 }
639}

References Nektar::SolverUtils::Coupling::m_couplingName, Nektar::SolverUtils::Coupling::m_evalField, m_sendHandle, Nektar::LibUtilities::Timer::Start(), Nektar::LibUtilities::Timer::Stop(), and Nektar::LibUtilities::Timer::TimePerTest().

Referenced by v_Send().

◆ SetupReceive()

void Nektar::SolverUtils::CouplingCwipi::SetupReceive ( )
private

Definition at line 236 of file CouplingCwipi.cpp.

237{
238 int oversamp = boost::lexical_cast<int>(m_config["OVERSAMPLE"]);
239
242 recvGraph->SetExpansionInfoToPointOrder(
243 oversamp + m_evalField->GetExp(0)->GetNumPoints(0));
244
245 // TODO: DeclareCoeffPhysArrays
247 m_evalField->GetSession(), recvGraph, "DefaultVar");
248
249 m_oldFields = Array<OneD, Array<OneD, NekDouble>>(m_nRecvVars);
250 m_newFields = Array<OneD, Array<OneD, NekDouble>>(m_nRecvVars);
251 for (int i = 0; i < m_nRecvVars; ++i)
252 {
253 m_oldFields[i] = Array<OneD, NekDouble>(m_evalField->GetTotPoints());
254 m_newFields[i] = Array<OneD, NekDouble>(m_evalField->GetTotPoints());
255 }
256
257 // define the quadrature points at which we want to receive data
258 m_nPoints = m_recvField->GetTotPoints();
259 Array<OneD, Array<OneD, NekDouble>> coords(3);
260 coords[0] = Array<OneD, NekDouble>(m_nPoints);
261 coords[1] = Array<OneD, NekDouble>(m_nPoints);
262 coords[2] = Array<OneD, NekDouble>(m_nPoints);
263 m_recvField->GetCoords(coords[0], coords[1], coords[2]);
264
265 m_points = (double *)malloc(sizeof(double) * 3 * m_nPoints);
266 ASSERTL1(m_points != nullptr, "malloc failed for m_points");
267
268 for (int i = 0; i < m_nPoints; ++i)
269 {
270 m_points[3 * i + 0] = double(coords[0][i]);
271
272 if (m_spacedim > 1)
273 {
274 m_points[3 * i + 1] = double(coords[1][i]);
275 }
276 else
277 {
278 m_points[3 * i + 1] = 0.0;
279 }
280
281 if (m_spacedim > 2)
282 {
283 m_points[3 * i + 2] = double(coords[2][i]);
284 }
285 else
286 {
287 m_points[3 * i + 2] = 0.0;
288 }
289 }
290
291 cwipi_set_points_to_locate(m_couplingName.c_str(), m_nPoints, m_points);
292
293 m_rValsInterl = (double *)malloc(sizeof(double) * m_nPoints * m_nRecvVars);
294 ASSERTL1(m_rValsInterl != nullptr, "malloc failed for m_rValsInterl");
295}
Array< OneD, Array< OneD, NekDouble > > m_oldFields
Definition: CouplingCwipi.h:96
Array< OneD, Array< OneD, NekDouble > > m_newFields
Definition: CouplingCwipi.h:97
static MeshGraphSharedPtr Read(const LibUtilities::SessionReaderSharedPtr pSession, LibUtilities::DomainRangeShPtr rng=LibUtilities::NullDomainRangeShPtr, bool fillGraph=true, SpatialDomains::MeshGraphSharedPtr partitionedGraph=nullptr)
Definition: MeshGraph.cpp:115

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), ASSERTL1, Nektar::SolverUtils::Coupling::m_config, Nektar::SolverUtils::Coupling::m_couplingName, Nektar::SolverUtils::Coupling::m_evalField, m_newFields, m_nPoints, Nektar::SolverUtils::Coupling::m_nRecvVars, m_oldFields, m_points, m_recvField, m_rValsInterl, m_spacedim, and Nektar::SpatialDomains::MeshGraph::Read().

Referenced by v_Init().

◆ SetupSend()

void Nektar::SolverUtils::CouplingCwipi::SetupSend ( )
private

Definition at line 297 of file CouplingCwipi.cpp.

298{
299 // this array is never used because of our send callback method
300 m_sValsInterl = (double *)malloc(
301 sizeof(double) * m_evalField->GetGraph()->GetNvertices() * m_nSendVars);
302 ASSERTL1(m_sValsInterl != nullptr, "malloc failed for m_sValsInterl");
303 for (int i = 0; i < m_evalField->GetGraph()->GetNvertices() * m_nSendVars;
304 ++i)
305 {
306 m_sValsInterl[i] = i;
307 }
308
309 // register this coupling as sender
310 std::stringstream sst;
311 sst << m_spacedim << "," << m_evalField->GetGraph()->GetNvertices() << ","
312 << m_nSendVars;
313 SendCallbackMap[sst.str()] =
314 std::bind(&CouplingCwipi::SendCallback, this, std::placeholders::_1,
315 std::placeholders::_2);
316 cwipi_set_interpolation_function(m_couplingName.c_str(),
318}
static SOLVER_UTILS_EXPORT void InterpCallback(const int entities_dim, const int n_local_vertex, const int n_local_element, const int n_local_polhyedra, const int n_distant_point, const double local_coordinates[], const int local_connectivity_index[], const int local_connectivity[], const int local_polyhedra_face_index[], const int local_polyhedra_cell_to_face_connectivity[], const int local_polyhedra_face_connectivity_index[], const int local_polyhedra_face_connectivity[], const double distant_points_coordinates[], const int distant_points_location[], const float distant_points_distance[], const int distant_points_barycentric_coordinates_index[], const double distant_points_barycentric_coordinates[], const int stride, const cwipi_solver_type_t solver_type, const void *local_field, void *distant_field)
SOLVER_UTILS_EXPORT void SendCallback(Array< OneD, Array< OneD, NekDouble > > &interpField, Array< OneD, Array< OneD, NekDouble > > &distCoords)

References ASSERTL1, InterpCallback(), Nektar::SolverUtils::Coupling::m_couplingName, Nektar::SolverUtils::Coupling::m_evalField, Nektar::SolverUtils::Coupling::m_nSendVars, m_spacedim, m_sValsInterl, SendCallback(), and Nektar::SolverUtils::SendCallbackMap.

Referenced by v_Init().

◆ SetupSendInterpolation()

void Nektar::SolverUtils::CouplingCwipi::SetupSendInterpolation ( )
private

Definition at line 371 of file CouplingCwipi.cpp.

372{
373 const double *distCoords =
374 cwipi_get_distant_coordinates(m_couplingName.c_str());
375 int nPts = cwipi_get_n_distant_points(m_couplingName.c_str());
376
377 Array<OneD, Array<OneD, NekDouble>> local(3);
378 for (int i = 0; i < 3; ++i)
379 {
380 local[i] = Array<OneD, NekDouble>(m_evalField->GetTotPoints(), 0.0);
381 }
382 m_evalField->GetCoords(local[0], local[1], local[2]);
385
386 Array<OneD, Array<OneD, NekDouble>> dist(3);
387 for (int i = 0; i < 3; ++i)
388 {
389 dist[i] = Array<OneD, NekDouble>(nPts);
390 for (int j = 0; j < nPts; ++j)
391 {
392 dist[i][j] = distCoords[3 * j + i];
393 }
394 }
397
399 if (boost::to_upper_copy(m_config["SENDMETHOD"]) == "SHEPARD")
400 {
401 method = LibUtilities::eShepard;
402 }
403 m_sendInterpolator = MemoryManager<FieldUtils::Interpolator<std::vector<
404 MultiRegions::ExpListSharedPtr>>>::AllocateSharedPtr(method);
405 m_sendInterpolator->CalcWeights(locatPts, distPts);
406 m_sendInterpolator->PrintStatistics();
407}
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::LibUtilities::eNearestNeighbour, Nektar::LibUtilities::eShepard, Nektar::SolverUtils::Coupling::m_config, Nektar::SolverUtils::Coupling::m_couplingName, Nektar::SolverUtils::Coupling::m_evalField, and m_sendInterpolator.

Referenced by SendCallback().

◆ v_Finalize()

void Nektar::SolverUtils::CouplingCwipi::v_Finalize ( void  )
overrideprotectedvirtual

Reimplemented from Nektar::SolverUtils::Coupling.

Definition at line 489 of file CouplingCwipi.cpp.

490{
491 cwipi_delete_coupling(m_couplingName.c_str());
492}

References Nektar::SolverUtils::Coupling::m_couplingName.

◆ v_Init()

void Nektar::SolverUtils::CouplingCwipi::v_Init ( )
overrideprotectedvirtual

Reimplemented from Nektar::SolverUtils::Coupling.

Definition at line 123 of file CouplingCwipi.cpp.

124{
126
127 ReadConfig(m_evalField->GetSession());
128
129 cwipi_add_local_int_control_parameter("nSendVars", m_nSendVars);
130 cwipi_add_local_int_control_parameter("nRecvVars", m_nRecvVars);
131 cwipi_add_local_string_control_parameter(
132 "recvFieldNames", m_config["RECEIVEVARIABLES"].c_str());
133 cwipi_add_local_string_control_parameter("sendFieldNames",
134 m_config["SENDVARIABLES"].c_str());
135
136 // MPI_TAG_UB is guaranteed to be at least 32767, so we make
137 // sure m_recvTag < 32767. Only caveat: m_recvTag is not guaranteed to be
138 // unique.
139 m_recvTag =
140 boost::hash<std::string>()(m_couplingName + m_config["REMOTENAME"] +
141 m_config["LOCALNAME"]) %
142 32767;
143
144 cwipi_add_local_int_control_parameter("receiveTag", m_recvTag);
145
146 m_spacedim = m_evalField->GetGraph()->GetSpaceDimension();
147
148 if (m_filtWidth > 0)
149 {
150 m_filtWidth = 2 * M_PI / m_filtWidth;
152 }
153
154 // Init Coupling
155 cwipi_solver_type_t solver_type = CWIPI_SOLVER_CELL_VERTEX;
156 NekDouble geom_tol = boost::lexical_cast<NekDouble>(m_config["GEOMTOL"]);
157 cwipi_create_coupling(
158 m_couplingName.c_str(), CWIPI_COUPLING_PARALLEL_WITH_PARTITIONING,
159 m_config["REMOTENAME"].c_str(), m_spacedim, geom_tol, CWIPI_STATIC_MESH,
160 solver_type, OUTPUT_FREQ, "Ensight Gold", "text");
161 cwipi_synchronize_control_parameter(m_config["REMOTENAME"].c_str());
162
163 if (m_evalField->GetComm()->GetRank() == 0 &&
164 m_evalField->GetSession()->DefinesCmdLineArgument("verbose"))
165 {
166 cwipi_dump_application_properties();
167 }
168
169 m_sendTag = cwipi_get_distant_int_control_parameter(
170 m_config["REMOTENAME"].c_str(), "receiveTag");
171
172 if (cwipi_has_int_parameter(m_config["REMOTENAME"].c_str(), "nRecvVars"))
173 {
174 int remoteNRecvVars = cwipi_get_distant_int_control_parameter(
175 m_config["REMOTENAME"].c_str(), "nRecvVars");
176 ASSERTL0(remoteNRecvVars == m_nSendVars,
177 "Number of local send vars different to remote received vars");
178 }
179
180 if (cwipi_has_int_parameter(m_config["REMOTENAME"].c_str(), "nSendVars"))
181 {
182 int remoteNSendVars = cwipi_get_distant_int_control_parameter(
183 m_config["REMOTENAME"].c_str(), "nSendVars");
184 ASSERTL0(remoteNSendVars == m_nRecvVars,
185 "Number of local receive vars different to remote sent vars");
186 }
187
188 AnnounceMesh();
189
190 if (m_nRecvVars > 0 && m_recvSteps > 0)
191 {
192 SetupReceive();
193 }
194
195 if (m_evalField->GetComm()->GetRank() == 0 &&
196 m_evalField->GetSession()->DefinesCmdLineArgument("verbose"))
197 {
198 cout << "locating..." << endl;
199 }
200 cwipi_locate(m_couplingName.c_str());
201
202 if (m_nSendVars > 0 && m_sendSteps > 0)
203 {
204 SetupSend();
205 }
206
207 if (m_nRecvVars > 0 && m_recvSteps > 0)
208 {
209 ReceiveStart();
210 }
211}
#define OUTPUT_FREQ
void ReadConfig(LibUtilities::SessionReaderSharedPtr session)
virtual SOLVER_UTILS_EXPORT void v_Init()
Definition: Coupling.cpp:59

References AnnounceMesh(), ASSERTL0, Nektar::SolverUtils::Coupling::m_config, Nektar::SolverUtils::Coupling::m_couplingName, Nektar::SolverUtils::Coupling::m_evalField, m_filtWidth, Nektar::SolverUtils::Coupling::m_nRecvVars, Nektar::SolverUtils::Coupling::m_nSendVars, Nektar::SolverUtils::Coupling::m_recvSteps, m_recvTag, Nektar::SolverUtils::Coupling::m_sendSteps, m_sendTag, m_spacedim, OUTPUT_FREQ, ReadConfig(), ReceiveStart(), SetupReceive(), SetupSend(), and Nektar::SolverUtils::Coupling::v_Init().

◆ v_Receive()

void Nektar::SolverUtils::CouplingCwipi::v_Receive ( const int  step,
const NekDouble  time,
Array< OneD, Array< OneD, NekDouble > > &  field,
std::vector< std::string > &  varNames 
)
overrideprotectedvirtual

Implements Nektar::SolverUtils::Coupling.

Definition at line 666 of file CouplingCwipi.cpp.

669{
670 if (m_nRecvVars < 1 || m_recvSteps < 1)
671 {
672 return;
673 }
674
675 Array<OneD, Array<OneD, NekDouble>> recvFields(m_nRecvVars);
676 vector<int> recvVarsToVars =
678 ASSERTL1(m_nRecvVars == recvVarsToVars.size(), "field size mismatch");
679 for (int i = 0; i < recvVarsToVars.size(); ++i)
680 {
681 recvFields[i] = field[recvVarsToVars[i]];
682 }
683
684 int nq = m_evalField->GetTotPoints();
685
686 // make sure we have sensible data in old/new recvFields the first time this
687 // method is called
688 if (m_lastReceive < 0)
689 {
690 for (int i = 0; i < m_nRecvVars; ++i)
691 {
692 Vmath::Vcopy(nq, recvFields[i], 1, m_oldFields[i], 1);
693 Vmath::Vcopy(nq, recvFields[i], 1, m_newFields[i], 1);
694 }
695 }
696
697 if (step >= m_lastReceive + m_recvSteps)
698 {
699 for (int i = 0; i < m_nRecvVars; ++i)
700 {
701 Vmath::Vcopy(nq, m_newFields[i], 1, m_oldFields[i], 1);
702 }
703
704 ReceiveCwipi(step, time, m_newFields);
705 }
706
707 NekDouble fact =
709 for (int i = 0; i < m_nRecvVars; ++i)
710 {
711 Vmath::Svtsvtp(nq, fact, m_newFields[i], 1, (1 - fact), m_oldFields[i],
712 1, recvFields[i], 1);
713 }
714}
void ReceiveCwipi(const int step, const NekDouble time, Array< OneD, Array< OneD, NekDouble > > &field)
SOLVER_UTILS_EXPORT std::vector< int > GenerateVariableMapping(std::vector< std::string > &vars, std::vector< std::string > &transVars)
Definition: Coupling.cpp:132
void Svtsvtp(int n, const T alpha, const T *x, int incx, const T beta, const T *y, int incy, T *z, int incz)
Svtsvtp (scalar times vector plus scalar times vector):
Definition: Vmath.hpp:473
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825

References ASSERTL1, Nektar::SolverUtils::Coupling::GenerateVariableMapping(), Nektar::SolverUtils::Coupling::m_evalField, m_lastReceive, m_newFields, Nektar::SolverUtils::Coupling::m_nRecvVars, m_oldFields, Nektar::SolverUtils::Coupling::m_recvFieldNames, Nektar::SolverUtils::Coupling::m_recvSteps, ReceiveCwipi(), Vmath::Svtsvtp(), and Vmath::Vcopy().

◆ v_Send()

void Nektar::SolverUtils::CouplingCwipi::v_Send ( const int  step,
const NekDouble  time,
const Array< OneD, const Array< OneD, NekDouble > > &  field,
std::vector< std::string > &  varNames 
)
overrideprotectedvirtual

Implements Nektar::SolverUtils::Coupling.

Definition at line 570 of file CouplingCwipi.cpp.

574{
575 if (m_nSendVars < 1 || m_sendSteps < 1)
576 {
577 return;
578 }
579
580 if (step >= m_lastSend + m_sendSteps)
581 {
582 SendComplete();
583
584 m_lastSend = step;
585
586 vector<int> sendVarsToVars =
588 m_sendField = Array<OneD, Array<OneD, NekDouble>>(m_nSendVars);
589 for (int i = 0; i < sendVarsToVars.size(); ++i)
590 {
591 m_sendField[i] = field[sendVarsToVars[i]];
592 }
593
594 if (m_evalField->GetComm()->GetRank() == 0 &&
595 m_evalField->GetSession()->DefinesCmdLineArgument("verbose"))
596 {
597 cout << "sending fields at i = " << step << ", t = " << time
598 << endl;
599 }
600
601 LibUtilities::Timer timer1;
602 timer1.Start();
603
604 char sendFN[10];
605 strcpy(sendFN, "dummyName");
606
607 cwipi_issend(m_couplingName.c_str(), "ex1", m_sendTag, m_nSendVars,
608 step, time, sendFN, m_sValsInterl, &m_sendHandle);
609 timer1.Stop();
610
611 if (m_evalField->GetComm()->GetRank() == 0 &&
612 m_evalField->GetSession()->DefinesCmdLineArgument("verbose"))
613 {
614 cout << "Send total time: " << timer1.TimePerTest(1) << endl;
615 }
616 }
617}
std::vector< std::string > m_sendFieldNames
Definition: Coupling.h:114

References Nektar::SolverUtils::Coupling::GenerateVariableMapping(), Nektar::SolverUtils::Coupling::m_couplingName, Nektar::SolverUtils::Coupling::m_evalField, m_lastSend, Nektar::SolverUtils::Coupling::m_nSendVars, m_sendField, Nektar::SolverUtils::Coupling::m_sendFieldNames, m_sendHandle, Nektar::SolverUtils::Coupling::m_sendSteps, m_sendTag, m_sValsInterl, SendComplete(), Nektar::LibUtilities::Timer::Start(), Nektar::LibUtilities::Timer::Stop(), and Nektar::LibUtilities::Timer::TimePerTest().

Member Data Documentation

◆ className

std::string Nektar::SolverUtils::CouplingCwipi::className
static
Initial value:
=
GetCouplingFactory().RegisterCreatorFunction("Cwipi", CouplingCwipi::create,
"Cwipi Coupling")
static CouplingSharedPtr create(MultiRegions::ExpListSharedPtr field)
Creates an instance of this class.
Definition: CouplingCwipi.h:56
CouplingFactory & GetCouplingFactory()
Declaration of the Coupling factory singleton.
Definition: Coupling.cpp:42

Definition at line 53 of file CouplingCwipi.h.

◆ m_connec

int* Nektar::SolverUtils::CouplingCwipi::m_connec
protected

Definition at line 113 of file CouplingCwipi.h.

Referenced by AddElementsToMesh(), AnnounceMesh(), and ~CouplingCwipi().

◆ m_connecIdx

int* Nektar::SolverUtils::CouplingCwipi::m_connecIdx
protected

Definition at line 112 of file CouplingCwipi.h.

Referenced by AddElementsToMesh(), AnnounceMesh(), and ~CouplingCwipi().

◆ m_coords

double* Nektar::SolverUtils::CouplingCwipi::m_coords
protected

Definition at line 111 of file CouplingCwipi.h.

Referenced by AddElementsToMesh(), AnnounceMesh(), and ~CouplingCwipi().

◆ m_extrapInterpolator

std::shared_ptr< FieldUtils::Interpolator<std::vector<MultiRegions::ExpListSharedPtr> > > Nektar::SolverUtils::CouplingCwipi::m_extrapInterpolator
protected

Definition at line 125 of file CouplingCwipi.h.

Referenced by ExtrapolateFields().

◆ m_filtWidth

NekDouble Nektar::SolverUtils::CouplingCwipi::m_filtWidth
protected

Definition at line 90 of file CouplingCwipi.h.

Referenced by ReadConfig(), ReceiveCwipi(), and v_Init().

◆ m_lastReceive

int Nektar::SolverUtils::CouplingCwipi::m_lastReceive
protected

Definition at line 106 of file CouplingCwipi.h.

Referenced by ReceiveCwipi(), and v_Receive().

◆ m_lastSend

int Nektar::SolverUtils::CouplingCwipi::m_lastSend
protected

Definition at line 105 of file CouplingCwipi.h.

Referenced by v_Send().

◆ m_newFields

Array<OneD, Array<OneD, NekDouble> > Nektar::SolverUtils::CouplingCwipi::m_newFields
protected

Definition at line 97 of file CouplingCwipi.h.

Referenced by SetupReceive(), and v_Receive().

◆ m_nPoints

int Nektar::SolverUtils::CouplingCwipi::m_nPoints
protected

Definition at line 109 of file CouplingCwipi.h.

Referenced by ExtrapolateFields(), ReceiveCwipi(), and SetupReceive().

◆ m_oldFields

Array<OneD, Array<OneD, NekDouble> > Nektar::SolverUtils::CouplingCwipi::m_oldFields
protected

Definition at line 96 of file CouplingCwipi.h.

Referenced by SetupReceive(), and v_Receive().

◆ m_points

double* Nektar::SolverUtils::CouplingCwipi::m_points
protected

Definition at line 110 of file CouplingCwipi.h.

Referenced by SetupReceive(), and ~CouplingCwipi().

◆ m_recvField

MultiRegions::ExpListSharedPtr Nektar::SolverUtils::CouplingCwipi::m_recvField
protected

Definition at line 94 of file CouplingCwipi.h.

Referenced by DumpRawFields(), ExtrapolateFields(), ReceiveCwipi(), and SetupReceive().

◆ m_recvHandle

int Nektar::SolverUtils::CouplingCwipi::m_recvHandle
protected

Definition at line 100 of file CouplingCwipi.h.

Referenced by ReceiveCwipi(), and ReceiveStart().

◆ m_recvTag

int Nektar::SolverUtils::CouplingCwipi::m_recvTag
protected

Definition at line 103 of file CouplingCwipi.h.

Referenced by ReceiveStart(), and v_Init().

◆ m_rValsInterl

double* Nektar::SolverUtils::CouplingCwipi::m_rValsInterl
protected

Definition at line 114 of file CouplingCwipi.h.

Referenced by ReceiveCwipi(), ReceiveStart(), SetupReceive(), and ~CouplingCwipi().

◆ m_sendField

Array<OneD, Array<OneD, NekDouble> > Nektar::SolverUtils::CouplingCwipi::m_sendField
protected

Definition at line 92 of file CouplingCwipi.h.

Referenced by EvaluateFields(), GetSendField(), SendCallback(), and v_Send().

◆ m_sendHandle

int Nektar::SolverUtils::CouplingCwipi::m_sendHandle
protected

Definition at line 99 of file CouplingCwipi.h.

Referenced by SendComplete(), and v_Send().

◆ m_sendInterpolator

std::shared_ptr< FieldUtils::Interpolator<std::vector<MultiRegions::ExpListSharedPtr> > > Nektar::SolverUtils::CouplingCwipi::m_sendInterpolator
protected

Definition at line 121 of file CouplingCwipi.h.

Referenced by SendCallback(), and SetupSendInterpolation().

◆ m_sendTag

int Nektar::SolverUtils::CouplingCwipi::m_sendTag
protected

Definition at line 102 of file CouplingCwipi.h.

Referenced by v_Init(), and v_Send().

◆ m_spacedim

int Nektar::SolverUtils::CouplingCwipi::m_spacedim
protected

◆ m_sValsInterl

double* Nektar::SolverUtils::CouplingCwipi::m_sValsInterl
protected

Definition at line 115 of file CouplingCwipi.h.

Referenced by SetupSend(), v_Send(), and ~CouplingCwipi().

◆ m_vertMap

std::map<int, int> Nektar::SolverUtils::CouplingCwipi::m_vertMap
protected

Definition at line 117 of file CouplingCwipi.h.

Referenced by AddElementsToMesh().