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)
 
virtual SOLVER_UTILS_EXPORT ~CouplingCwipi ()
 
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

virtual SOLVER_UTILS_EXPORT void v_Init () override
 
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) override
 
virtual SOLVER_UTILS_EXPORT void v_Receive (const int step, const NekDouble time, Array< OneD, Array< OneD, NekDouble >> &field, std::vector< std::string > &varNames) override
 
virtual 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)
 
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 51 of file CouplingCwipi.h.

Constructor & Destructor Documentation

◆ CouplingCwipi()

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

Definition at line 114 of file CouplingCwipi.cpp.

115  : Coupling(field), m_sendHandle(-1), m_recvHandle(-1), m_lastSend(-1E6),
116  m_lastReceive(-1E6), m_points(NULL), m_coords(NULL), m_connecIdx(NULL),
117  m_connec(NULL), m_rValsInterl(NULL), m_sValsInterl(NULL)
118 {
119  // defaults
120  m_config["GEOMTOL"] = "0.1";
121  m_config["LOCALNAME"] = "nektar";
122  m_config["REMOTENAME"] = "precise";
123  m_config["OVERSAMPLE"] = "0";
124  m_config["FILTERWIDTH"] = "-1";
125  m_config["DUMPRAW"] = "0";
126  m_config["SENDMETHOD"] = "EVALUATE";
127  m_config["NOTLOCMETHOD"] = "KEEP";
128 }
CouplingConfigMap m_config
Definition: Coupling.h:111
SOLVER_UTILS_EXPORT Coupling(MultiRegions::ExpListSharedPtr field)
Definition: Coupling.cpp:50

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

◆ ~CouplingCwipi()

Nektar::SolverUtils::CouplingCwipi::~CouplingCwipi ( )
virtual

Definition at line 220 of file CouplingCwipi.cpp.

221 {
222  free(m_coords);
223  free(m_points);
224  free(m_connec);
225  free(m_connecIdx);
226  free(m_rValsInterl);
227  free(m_sValsInterl);
228 }

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 502 of file CouplingCwipi.cpp.

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

References m_connec, m_connecIdx, m_coords, and m_vertMap.

Referenced by AnnounceMesh().

◆ AnnounceMesh()

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

Definition at line 416 of file CouplingCwipi.cpp.

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

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 58 of file CouplingCwipi.h.

59  {
62  p->Init();
63  return p;
64  }
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 1008 of file CouplingCwipi.cpp.

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

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 327 of file CouplingCwipi.cpp.

330 {
331  int nOutPts = distCoords.size();
332 
333  Array<OneD, NekDouble> Lcoords(m_spacedim, 0.0);
334  for (int i = 0; i < nOutPts; ++i)
335  {
336  // Obtain Element and LocalCoordinate to interpolate
337  int elmtid = -1;
339  while (elmtid < 0 && tol <= 1E3 * NekConstants::kNekZeroTol)
340  {
341  elmtid = m_evalField->GetExpIndex(distCoords[i], Lcoords, tol);
342  tol *= 2;
343  }
344  if (tol > 2 * NekConstants::kNekZeroTol)
345  {
346  for (int j = 0; j < m_spacedim; ++j)
347  {
348  if (Lcoords[j] < -1 - 0.75 * NekConstants::kNekZeroTol)
349  {
350  Lcoords[j] = -1;
351  }
352  if (Lcoords[j] > 1 + 0.75 * NekConstants::kNekZeroTol)
353  {
354  Lcoords[j] = 1;
355  }
356  }
357  }
358 
359  ASSERTL0(elmtid >= 0,
360  "no element found for (" +
361  boost::lexical_cast<string>(distCoords[i][0]) + ", " +
362  boost::lexical_cast<string>(distCoords[i][1]) + ", " +
363  boost::lexical_cast<string>(distCoords[i][2]) + ")");
364 
365  int offset = m_evalField->GetPhys_Offset(elmtid);
366 
367  for (int f = 0; f < m_nSendVars; ++f)
368  {
369  NekDouble value = m_evalField->GetExp(elmtid)->StdPhysEvaluate(
370  Lcoords, m_sendField[f] + offset);
371 
372  ASSERTL0(!(boost::math::isnan)(value), "new value is not a number");
373  interpField[f][i] = value;
374  }
375  }
376 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
Array< OneD, Array< OneD, NekDouble > > m_sendField
Definition: CouplingCwipi.h:94
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 882 of file CouplingCwipi.cpp.

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

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 143 of file CouplingCwipi.h.

144  {
145  return m_sendField[i][j];
146  }

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 62 of file CouplingCwipi.cpp.

76 {
77  boost::ignore_unused(
78  n_local_element, n_local_polyhedra, local_coordinates,
79  local_connectivity_index, local_connectivity,
80  local_polyhedra_face_index, local_polyhedra_cell_to_face_connectivity,
81  local_polyhedra_face_connectivity_index,
82  local_polyhedra_face_connectivity, distant_points_location,
83  distant_points_distance, distant_points_barycentric_coordinates_index,
84  distant_points_barycentric_coordinates, solver_type, local_field);
85 
86  Array<OneD, Array<OneD, NekDouble>> interpField(stride);
87 
88  Array<OneD, Array<OneD, NekDouble>> distCoords(n_distant_point);
89  for (int i = 0; i < n_distant_point; ++i)
90  {
91  distCoords[i] = Array<OneD, NekDouble>(3);
92  for (int j = 0; j < 3; ++j)
93  {
94  distCoords[i][j] = distant_points_coordinates[3 * i + j];
95  }
96  }
97 
98  std::stringstream sst;
99  sst << entities_dim << "," << n_local_vertex << "," << stride;
100  SendCallbackMap[sst.str()](interpField, distCoords);
101 
102  ASSERTL0(interpField.size() == stride, "size mismatch");
103  ASSERTL0(interpField[0].size() == n_distant_point, "size mismatch");
104 
105  for (int i = 0; i < n_distant_point; i++)
106  {
107  for (int j = 0; j < stride; ++j)
108  {
109  ((double *)distant_field)[i * stride + j] = interpField[j][i];
110  }
111  }
112 }
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 230 of file CouplingCwipi.cpp.

231 {
232  ASSERTL0(session->DefinesElement("Nektar/Coupling"),
233  "No Coupling config found");
234 
235  m_config["LOCALNAME"] = session->GetCmdLineArgument<std::string>("cwipi");
236 
237  TiXmlElement *vCoupling = session->GetElement("Nektar/Coupling");
238  ASSERTL0(vCoupling, "Invalid Coupling config");
239 
240  m_filtWidth = boost::lexical_cast<NekDouble>(m_config["FILTERWIDTH"]);
241 }

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 723 of file CouplingCwipi.cpp.

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

References ASSERTL1, DumpRawFields(), ExtrapolateFields(), 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, 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 648 of file CouplingCwipi.cpp.

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

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 544 of file CouplingCwipi.cpp.

547 {
548  ASSERTL0(interpField.size() == m_nSendVars, "size mismatch");
549 
550  for (int i = 0; i < m_nSendVars; ++i)
551  {
552  interpField[i] = Array<OneD, NekDouble>(distCoords.size());
553  }
554 
555  if (boost::to_upper_copy(m_config["SENDMETHOD"]) == "NEARESTNEIGHBOUR" ||
556  boost::to_upper_copy(m_config["SENDMETHOD"]) == "SHEPARD")
557  {
558  if (!m_sendInterpolator)
559  {
561  }
562 
565  0, m_sendField);
568  0, interpField);
569  m_sendInterpolator->Interpolate(ptsIn, ptsOut);
570  }
571  else
572  {
573  EvaluateFields(interpField, distCoords);
574  }
575 }
void EvaluateFields(Array< OneD, Array< OneD, NekDouble >> interpField, Array< OneD, Array< OneD, NekDouble >> distCoords)
std::shared_ptr< FieldUtils::Interpolator< std::vector< MultiRegions::ExpListSharedPtr > > > m_sendInterpolator

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 626 of file CouplingCwipi.cpp.

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

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 243 of file CouplingCwipi.cpp.

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

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 304 of file CouplingCwipi.cpp.

305 {
306  // this array is never used because of our send callback method
307  m_sValsInterl = (double *)malloc(
308  sizeof(double) * m_evalField->GetGraph()->GetNvertices() * m_nSendVars);
309  ASSERTL1(m_sValsInterl != NULL, "malloc failed for m_sValsInterl");
310  for (int i = 0; i < m_evalField->GetGraph()->GetNvertices() * m_nSendVars;
311  ++i)
312  {
313  m_sValsInterl[i] = i;
314  }
315 
316  // register this coupling as sender
317  std::stringstream sst;
318  sst << m_spacedim << "," << m_evalField->GetGraph()->GetNvertices() << ","
319  << m_nSendVars;
320  SendCallbackMap[sst.str()] =
321  std::bind(&CouplingCwipi::SendCallback, this, std::placeholders::_1,
322  std::placeholders::_2);
323  cwipi_set_interpolation_function(m_couplingName.c_str(),
325 }
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 378 of file CouplingCwipi.cpp.

379 {
380  const double *distCoords =
381  cwipi_get_distant_coordinates(m_couplingName.c_str());
382  int nPts = cwipi_get_n_distant_points(m_couplingName.c_str());
383 
384  Array<OneD, Array<OneD, NekDouble>> local(3);
385  for (int i = 0; i < 3; ++i)
386  {
387  local[i] = Array<OneD, NekDouble>(m_evalField->GetTotPoints(), 0.0);
388  }
389  m_evalField->GetCoords(local[0], local[1], local[2]);
392 
393  Array<OneD, Array<OneD, NekDouble>> dist(3);
394  for (int i = 0; i < 3; ++i)
395  {
396  dist[i] = Array<OneD, NekDouble>(nPts);
397  for (int j = 0; j < nPts; ++j)
398  {
399  dist[i][j] = distCoords[3 * j + i];
400  }
401  }
404 
406  if (boost::to_upper_copy(m_config["SENDMETHOD"]) == "SHEPARD")
407  {
408  method = LibUtilities::eShepard;
409  }
410  m_sendInterpolator = MemoryManager<FieldUtils::Interpolator<std::vector<
411  MultiRegions::ExpListSharedPtr>>>::AllocateSharedPtr(method);
412  m_sendInterpolator->CalcWeights(locatPts, distPts);
413  m_sendInterpolator->PrintStatistics();
414 }
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 496 of file CouplingCwipi.cpp.

497 {
498  cwipi_delete_coupling(m_couplingName.c_str());
499 }

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

◆ v_Init()

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

Reimplemented from Nektar::SolverUtils::Coupling.

Definition at line 130 of file CouplingCwipi.cpp.

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

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 673 of file CouplingCwipi.cpp.

676 {
677  if (m_nRecvVars < 1 || m_recvSteps < 1)
678  {
679  return;
680  }
681 
682  Array<OneD, Array<OneD, NekDouble>> recvFields(m_nRecvVars);
683  vector<int> recvVarsToVars =
685  ASSERTL1(m_nRecvVars == recvVarsToVars.size(), "field size mismatch");
686  for (int i = 0; i < recvVarsToVars.size(); ++i)
687  {
688  recvFields[i] = field[recvVarsToVars[i]];
689  }
690 
691  int nq = m_evalField->GetTotPoints();
692 
693  // make sure we have sensible data in old/new recvFields the first time this
694  // method is called
695  if (m_lastReceive < 0)
696  {
697  for (int i = 0; i < m_nRecvVars; ++i)
698  {
699  Vmath::Vcopy(nq, recvFields[i], 1, m_oldFields[i], 1);
700  Vmath::Vcopy(nq, recvFields[i], 1, m_newFields[i], 1);
701  }
702  }
703 
704  if (step >= m_lastReceive + m_recvSteps)
705  {
706  for (int i = 0; i < m_nRecvVars; ++i)
707  {
708  Vmath::Vcopy(nq, m_newFields[i], 1, m_oldFields[i], 1);
709  }
710 
711  ReceiveCwipi(step, time, m_newFields);
712  }
713 
714  NekDouble fact =
716  for (int i = 0; i < m_nRecvVars; ++i)
717  {
718  Vmath::Svtsvtp(nq, fact, m_newFields[i], 1, (1 - fact), m_oldFields[i],
719  1, recvFields[i], 1);
720  }
721 }
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:134
void Svtsvtp(int n, const T alpha, const T *x, int incx, const T beta, const T *y, int incy, T *z, int incz)
svtvvtp (scalar times vector plus scalar times vector):
Definition: Vmath.cpp:751
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1255

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 577 of file CouplingCwipi.cpp.

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

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:58
CouplingFactory & GetCouplingFactory()
Declaration of the Coupling factory singleton.
Definition: Coupling.cpp:44

Definition at line 55 of file CouplingCwipi.h.

◆ m_connec

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

Definition at line 115 of file CouplingCwipi.h.

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

◆ m_connecIdx

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

Definition at line 114 of file CouplingCwipi.h.

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

◆ m_coords

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

Definition at line 113 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 127 of file CouplingCwipi.h.

Referenced by ExtrapolateFields().

◆ m_filtWidth

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

Definition at line 92 of file CouplingCwipi.h.

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

◆ m_lastReceive

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

Definition at line 108 of file CouplingCwipi.h.

Referenced by ReceiveCwipi(), and v_Receive().

◆ m_lastSend

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

Definition at line 107 of file CouplingCwipi.h.

Referenced by v_Send().

◆ m_newFields

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

Definition at line 99 of file CouplingCwipi.h.

Referenced by SetupReceive(), and v_Receive().

◆ m_nPoints

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

Definition at line 111 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 98 of file CouplingCwipi.h.

Referenced by SetupReceive(), and v_Receive().

◆ m_points

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

Definition at line 112 of file CouplingCwipi.h.

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

◆ m_recvField

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

Definition at line 96 of file CouplingCwipi.h.

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

◆ m_recvHandle

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

Definition at line 102 of file CouplingCwipi.h.

Referenced by ReceiveCwipi(), and ReceiveStart().

◆ m_recvTag

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

Definition at line 105 of file CouplingCwipi.h.

Referenced by ReceiveStart(), and v_Init().

◆ m_rValsInterl

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

Definition at line 116 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 94 of file CouplingCwipi.h.

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

◆ m_sendHandle

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

Definition at line 101 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 123 of file CouplingCwipi.h.

Referenced by SendCallback(), and SetupSendInterpolation().

◆ m_sendTag

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

Definition at line 104 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 117 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 119 of file CouplingCwipi.h.

Referenced by AddElementsToMesh().