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 ()
 
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)
 
virtual SOLVER_UTILS_EXPORT void v_Receive (const int step, const NekDouble time, Array< OneD, Array< OneD, NekDouble > > &field, std::vector< std::string > &varNames)
 
virtual SOLVER_UTILS_EXPORT void v_Finalize ()
 
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
 
FieldUtils::InterpolatorSharedPtr m_sendInterpolator
 
FieldUtils::InterpolatorSharedPtr 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 126 of file CouplingCwipi.cpp.

127  : Coupling(field), m_sendHandle(-1), m_recvHandle(-1), m_lastSend(-1E6),
128  m_lastReceive(-1E6), m_points(NULL), m_coords(NULL), m_connecIdx(NULL),
129  m_connec(NULL), m_rValsInterl(NULL), m_sValsInterl(NULL)
130 {
131  // defaults
132  m_config["GEOMTOL"] = "0.1";
133  m_config["LOCALNAME"] = "nektar";
134  m_config["REMOTENAME"] = "precise";
135  m_config["OVERSAMPLE"] = "0";
136  m_config["FILTERWIDTH"] = "-1";
137  m_config["DUMPRAW"] = "0";
138  m_config["SENDMETHOD"] = "EVALUATE";
139  m_config["NOTLOCMETHOD"] = "KEEP";
140 }
CouplingConfigMap m_config
Definition: Coupling.h:113
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 236 of file CouplingCwipi.cpp.

237 {
238  free(m_coords);
239  free(m_points);
240  free(m_connec);
241  free(m_connecIdx);
242  free(m_rValsInterl);
243  free(m_sValsInterl);
244 }

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

524 {
525  // helper variables
526  Array<OneD, NekDouble> x(3);
528  int vertID;
529 
530  int kNverts = T::mapped_type::element_type::kNverts;
531 
532  // iterate over all elements
533  for (auto it = geom.begin(); it != geom.end(); it++)
534  {
535  // iterate over the elements vertices
536  for (int j = 0; j < kNverts; ++j)
537  {
538  vert = it->second->GetVertex(j);
539  vertID = vert->GetVid();
540 
541  // check if we already stored the vertex
542  if (m_vertMap.count(vertID) == 0)
543  {
544  // store the vertex
545  vert->GetCoords(x[0], x[1], x[2]);
546  m_coords[3 * coordsPos + 0] = double(x[0]);
547  m_coords[3 * coordsPos + 1] = double(x[1]);
548  m_coords[3 * coordsPos + 2] = double(x[2]);
549 
550  // store the vertex position in the m_coords array
551  m_vertMap[vertID] = coordsPos;
552  coordsPos++;
553  }
554 
555  m_connec[connecPos] = m_vertMap[vertID] + 1;
556  connecPos++;
557  }
558 
559  m_connecIdx[conidxPos + 1] = m_connecIdx[conidxPos] + kNverts;
560  conidxPos++;
561  }
562 }
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 434 of file CouplingCwipi.cpp.

435 {
437 
438  // get Elements
446  if (m_spacedim == 1)
447  {
448  seggeom = graph->GetAllSegGeoms();
449  }
450  else if (m_spacedim == 2)
451  {
452  trigeom = graph->GetAllTriGeoms();
453  quadgeom = graph->GetAllQuadGeoms();
454  }
455  else if (m_spacedim == 3)
456  {
457  tetgeom = graph->GetAllTetGeoms();
458  pyrgeom = graph->GetAllPyrGeoms();
459  prismgeom = graph->GetAllPrismGeoms();
460  hexgeom = graph->GetAllHexGeoms();
461  };
462 
463  int nVerts = graph->GetNvertices();
464  int nElts = seggeom.size() + trigeom.size() + quadgeom.size() +
465  tetgeom.size() + pyrgeom.size() + prismgeom.size() +
466  hexgeom.size();
467 
468  // allocate CWIPI arrays
469  m_coords = (double *)malloc(sizeof(double) * 3 * nVerts);
470  ASSERTL1(m_coords != NULL, "malloc failed for m_coords");
471  int tmp = 2 * seggeom.size() + 3 * trigeom.size() + 4 * quadgeom.size() +
472  4 * tetgeom.size() + 5 * pyrgeom.size() + 6 * prismgeom.size() +
473  8 * hexgeom.size();
474  m_connec = (int *)malloc(sizeof(int) * tmp);
475  ASSERTL1(m_connec != NULL, "malloc failed for m_connec");
476  m_connecIdx = (int *)malloc(sizeof(int) * (nElts + 1));
477  ASSERTL1(m_connecIdx != NULL, "malloc failed for m_connecIdx");
478 
479  m_connecIdx[0] = 0;
480  int coordsPos = 0;
481  int connecPos = 0;
482  int conidxPos = 0;
483 
484  AddElementsToMesh(seggeom, coordsPos, connecPos, conidxPos);
485  AddElementsToMesh(trigeom, coordsPos, connecPos, conidxPos);
486  AddElementsToMesh(quadgeom, coordsPos, connecPos, conidxPos);
487  AddElementsToMesh(tetgeom, coordsPos, connecPos, conidxPos);
488  AddElementsToMesh(pyrgeom, coordsPos, connecPos, conidxPos);
489  AddElementsToMesh(prismgeom, coordsPos, connecPos, conidxPos);
490  AddElementsToMesh(hexgeom, coordsPos, connecPos, conidxPos);
491 
492  // output the mesh in tecplot format. If this works, CWIPI will be able
493  // to process it, too
494  /*
495  cout << "VARIABLES = \"X\", \"Y\", \"Z\", \"U\"" << endl;
496  cout << "ZONE NODES=" << nVerts << ",ELEMENTS=" << nElts;
497  cout << ",DATAPACKING=POINT, ZONETYPE=FETETRAHEDRON" << endl;
498  for (int i = 0; i < nVerts; ++i)
499  {
500  cout << m_coords[3*i + 0] << " " << m_coords[3*i + 1] << " ";
501  cout << m_coords[3*i + 2] << " " << 1.0 << endl;
502  }
503  for (int i = 0; i < nElts; ++i)
504  {
505  cout << m_connec[i*4 + 0] << " " << m_connec[i*4 + 1] << " ";
506  cout << m_connec[i*4 + 2] << " " << m_connec[i*4 + 3] << endl;
507  }
508  */
509 
510  cwipi_define_mesh(
511  m_couplingName.c_str(), nVerts, nElts, m_coords, m_connecIdx, m_connec);
512 }
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:250
void AddElementsToMesh(T geom, int &coordsPos, int &connecPos, int &conidxPos)
MultiRegions::ExpListSharedPtr m_evalField
Definition: Coupling.h:115
std::map< int, TriGeomSharedPtr > TriGeomMap
Definition: TriGeom.h:59
std::map< int, PyrGeomSharedPtr > PyrGeomMap
Definition: PyrGeom.h:75
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: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 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 1052 of file CouplingCwipi.cpp.

1054 {
1055  LibUtilities::Timer timer1;
1056  timer1.Start();
1057 
1058 #if (defined _WIN32 && _MSC_VER < 1900)
1059  // We need this to make sure boost::format has always
1060  // two digits in the exponents of Scientific notation.
1061  unsigned int old_exponent_format;
1062  old_exponent_format = _set_output_format(_TWO_DIGIT_EXPONENT);
1063  filename = boost::str(boost::format(filename) % m_time);
1064  _set_output_format(old_exponent_format);
1065 #else
1066  std::string filename =
1067  boost::str(boost::format(m_config["DUMPRAW"]) % time);
1068 #endif
1069 
1070  Array<OneD, Array<OneD, NekDouble> > tmp(m_nRecvVars + m_spacedim);
1071  for (int i = 0; i < 3; ++i)
1072  {
1073  tmp[i] = Array<OneD, NekDouble>(m_recvField->GetTotPoints(), 0.0);
1074  }
1075  m_recvField->GetCoords(tmp[0], tmp[1], tmp[2]);
1076 
1077  for (int i = 0; i < m_nRecvVars; ++i)
1078  {
1079  tmp[m_spacedim + i] = rVals[i];
1080  }
1081 
1082  LibUtilities::CsvIO csvIO(m_evalField->GetSession()->GetComm());
1086  csvIO.Write(filename, rvPts);
1087 
1088  timer1.Stop();
1089  if (m_evalField->GetComm()->GetRank() == 0 &&
1090  m_evalField->GetSession()->DefinesCmdLineArgument("verbose"))
1091  {
1092  cout << "DumpRawFields total time: " << timer1.TimePerTest(1) << endl;
1093  }
1094 }
MultiRegions::ExpListSharedPtr m_recvField
std::vector< std::string > m_recvFieldNames
Definition: Coupling.h:122
std::shared_ptr< PtsField > PtsFieldSharedPtr
Definition: PtsField.h:182

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

348 {
349  int nOutPts = distCoords.size();
350 
351  Array<OneD, NekDouble> Lcoords(m_spacedim, 0.0);
352  for (int i = 0; i < nOutPts; ++i)
353  {
354  // Obtain Element and LocalCoordinate to interpolate
355  int elmtid = -1;
357  while (elmtid < 0 && tol <= 1E3 * NekConstants::kNekZeroTol)
358  {
359  elmtid = m_evalField->GetExpIndex(distCoords[i], Lcoords, tol);
360  tol *= 2;
361  }
362  if (tol > 2 * NekConstants::kNekZeroTol)
363  {
364  for (int j = 0; j < m_spacedim; ++j)
365  {
366  if (Lcoords[j] < -1 - 0.75 * NekConstants::kNekZeroTol)
367  {
368  Lcoords[j] = -1;
369  }
370  if (Lcoords[j] > 1 + 0.75 * NekConstants::kNekZeroTol)
371  {
372  Lcoords[j] = 1;
373  }
374  }
375  }
376 
377  ASSERTL0(elmtid >= 0,
378  "no element found for (" +
379  boost::lexical_cast<string>(distCoords[i][0]) + ", " +
380  boost::lexical_cast<string>(distCoords[i][1]) + ", " +
381  boost::lexical_cast<string>(distCoords[i][2]) + ")");
382 
383  int offset = m_evalField->GetPhys_Offset(elmtid);
384 
385  for (int f = 0; f < m_nSendVars; ++f)
386  {
387  NekDouble value = m_evalField->GetExp(elmtid)->StdPhysEvaluate(
388  Lcoords, m_sendField[f] + offset);
389 
390  ASSERTL0(!(boost::math::isnan)(value), "new value is not a number");
391  interpField[f][i] = value;
392  }
393  }
394 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
Array< OneD, Array< OneD, NekDouble > > m_sendField
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 926 of file CouplingCwipi.cpp.

928 {
929  LibUtilities::Timer timer1, timer2;
930  timer1.Start();
931 
932  int totvars = 3 + m_nRecvVars;
933  int nNotLoc = notLoc.size();
934  int nranks = m_evalField->GetSession()->GetComm()->GetSize();
935 
936  Array<OneD, int> thisNLoc(1, m_nPoints - nNotLoc);
937  Array<OneD, int> sizeMap(nranks);
938  m_evalField->GetSession()->GetComm()->AllGather(thisNLoc, sizeMap);
939 
940  Array<OneD, int> offsetMap(nranks);
941  offsetMap[0] = 0;
942  for (int i = 1; i < nranks; ++i)
943  {
944  offsetMap[i] = offsetMap[i - 1] + sizeMap[i - 1];
945  }
946  int totNLoc = offsetMap[nranks - 1] + sizeMap[nranks - 1];
947 
948  Array<OneD, Array<OneD, NekDouble> > allVals(totvars);
949  for (int i = 0; i < 3; ++i)
950  {
951  allVals[i] = Array<OneD, NekDouble>(m_nPoints);
952  }
953  m_recvField->GetCoords(allVals[0], allVals[1], allVals[2]);
954 
955  if (m_spacedim < 3)
956  {
957  Vmath::Zero(m_nPoints, allVals[2], 1);
958  }
959  if (m_spacedim < 2)
960  {
961  Vmath::Zero(m_nPoints, allVals[1], 1);
962  }
963 
964  for (int i = 0; i < m_nRecvVars; ++i)
965  {
966  allVals[3 + i] = rVals[i];
967  }
968 
969  Array<OneD, Array<OneD, NekDouble> > locatedVals(totvars);
970  for (int i = 0; i < totvars; ++i)
971  {
972  locatedVals[i] = Array<OneD, NekDouble>(totNLoc, -42.0);
973  }
974 
975  // only copy points from allVals to locatedVals that were located
976  int offset = offsetMap[m_evalField->GetSession()->GetComm()->GetRank()];
977  for (int i = 0, intPos = 0, locPos = 0; i < m_nPoints; ++i)
978  {
979  if (locPos < nNotLoc && notLoc[locPos] == i)
980  {
981  // do nothing
982  locPos++;
983  }
984  else
985  {
986  for (int k = 0; k < totvars; ++k)
987  {
988  locatedVals[k][offset + intPos] = allVals[k][i];
989  }
990  intPos++;
991  }
992  }
993 
994  // send all located points to all ranks. This is probably horribly expensive
995  timer2.Start();
996  for (int i = 0; i < totvars; ++i)
997  {
998  m_evalField->GetSession()->GetComm()->AllGatherv(
999  locatedVals[i], sizeMap, offsetMap);
1000  }
1001  timer2.Stop();
1002 
1003  if (nNotLoc > 0)
1004  {
1005  LibUtilities::PtsFieldSharedPtr locatedPts =
1007  3, locatedVals);
1008 
1009  Array<OneD, Array<OneD, NekDouble> > notLocVals(totvars);
1010  for (int j = 0; j < totvars; ++j)
1011  {
1012  notLocVals[j] = Array<OneD, NekDouble>(nNotLoc);
1013  for (int i = 0; i < nNotLoc; ++i)
1014  {
1015  notLocVals[j][i] = allVals[j][notLoc[i]];
1016  }
1017  }
1020  3, notLocVals);
1021 
1022  // perform a nearest neighbour interpolation from locatedVals to the not
1023  // located rVals
1024  if (!m_extrapInterpolator)
1025  {
1029  m_extrapInterpolator->CalcWeights(locatedPts, notlocPts);
1030  m_extrapInterpolator->PrintStatistics();
1031  }
1032  m_extrapInterpolator->Interpolate(locatedPts, notlocPts);
1033 
1034  for (int j = 3; j < totvars; ++j)
1035  {
1036  for (int i = 0; i < nNotLoc; ++i)
1037  {
1038  allVals[j][notLoc[i]] = notLocVals[j][i];
1039  }
1040  }
1041  }
1042 
1043  timer1.Stop();
1044  if (m_evalField->GetComm()->GetRank() == 0 &&
1045  m_evalField->GetSession()->DefinesCmdLineArgument("verbose"))
1046  {
1047  cout << "ExtrapolateFields total time: " << timer1.TimePerTest(1);
1048  cout << " (AllGathervI: " << timer2.TimePerTest(1) << ")" << endl;
1049  }
1050 }
FieldUtils::InterpolatorSharedPtr m_extrapInterpolator
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:436

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

149  {
150  return m_sendField[i][j];
151  }

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.

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

247 {
248  ASSERTL0(session->DefinesElement("Nektar/Coupling"),
249  "No Coupling config found");
250 
251  m_config["LOCALNAME"] = session->GetCmdLineArgument<std::string>("cwipi");
252 
253  TiXmlElement *vCoupling = session->GetElement("Nektar/Coupling");
254  ASSERTL0(vCoupling, "Invalid Coupling config");
255 
256  m_filtWidth = boost::lexical_cast<NekDouble>(m_config["FILTERWIDTH"]);
257 }

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

769 {
770  ASSERTL1(m_nRecvVars == field.size(), "field size mismatch");
771 
772  if (m_nRecvVars < 1 || m_recvSteps < 1)
773  {
774  return;
775  }
776 
777  if (step >= m_lastReceive + m_recvSteps)
778  {
779  m_lastReceive = step;
780 
781  if (m_evalField->GetComm()->GetRank() == 0 &&
782  m_evalField->GetSession()->DefinesCmdLineArgument("verbose"))
783  {
784  cout << "waiting for receive at i = " << step << ", t = " << time
785  << endl;
786  }
787 
788  LibUtilities::Timer timer1, timer2, timer3;
789  timer1.Start();
790 
791  Array<OneD, NekDouble> tmpC(m_recvField->GetNcoeffs());
792  Array<OneD, Array<OneD, NekDouble> > rVals(m_nRecvVars);
793  for (int i = 0; i < m_nRecvVars; ++i)
794  {
795  rVals[i] = Array<OneD, NekDouble>(m_recvField->GetTotPoints());
796  }
797 
798  timer2.Start();
799  cwipi_wait_irecv(m_couplingName.c_str(), m_recvHandle);
800  timer2.Stop();
801 
802  // set to -1 so we know we can start receiving again
803  m_recvHandle = -1;
804 
805  int nNotLoc = cwipi_get_n_not_located_points(m_couplingName.c_str());
806  Array<OneD, int> notLoc;
807  if (nNotLoc != 0)
808  {
809  cout << "WARNING: relocating " << nNotLoc << " of " << m_nPoints
810  << " points" << endl;
811 
812  const int *tmp =
813  cwipi_get_not_located_points(m_couplingName.c_str());
814  notLoc = Array<OneD, int>(nNotLoc);
815  for (int i = 0; i < nNotLoc; ++i)
816  {
817  notLoc[i] = tmp[i] - 1;
818  }
819 
820  if (boost::to_upper_copy(m_config["NOTLOCMETHOD"]) == "KEEP")
821  {
822  // interpolate from m_evalField to m_recvField
823  for (int i = 0; i < m_nRecvVars; ++i)
824  {
825  m_evalField->FwdTrans(field[i], tmpC);
826  m_recvField->BwdTrans(tmpC, rVals[i]);
827  }
828  }
829  }
830 
831  for (int i = 0, locPos = 0, intPos = 0; i < m_nPoints; ++i)
832  {
833  if (locPos < nNotLoc && notLoc[locPos] == i)
834  {
835  // keep the original value of field[j][i]
836  locPos++;
837  }
838  else
839  {
840  for (int j = 0; j < m_nRecvVars; ++j)
841  {
842  rVals[j][i] = m_rValsInterl[intPos * m_nRecvVars + j];
843  }
844  intPos++;
845  }
846  }
847 
848  if (boost::to_upper_copy(m_config["NOTLOCMETHOD"]) == "EXTRAPOLATE")
849  {
850  int doExtrapolate = 0;
851  if (nNotLoc != 0)
852  {
853  doExtrapolate = 1;
854  }
855  m_evalField->GetSession()->GetComm()->AllReduce(
856  doExtrapolate, LibUtilities::ReduceMax);
857  if (doExtrapolate > 0)
858  {
859  ExtrapolateFields(rVals, notLoc);
860  }
861  }
862 
863  if (m_config["DUMPRAW"] != "0")
864  {
865  DumpRawFields(time, rVals);
866  }
867 
868  if (m_filtWidth > 0)
869  {
870  for (int i = 0; i < m_nRecvVars; ++i)
871  {
872  timer3.Start();
873 
874  Array<OneD, NekDouble> forcing(m_nPoints);
875 
876  Array<OneD, Array<OneD, NekDouble> > Velocity(m_spacedim);
877  for (int j = 0; j < m_spacedim; ++j)
878  {
879  Velocity[j] = Array<OneD, NekDouble>(m_nPoints, 0.0);
880  }
881 
882  Vmath::Smul(m_nPoints, -m_filtWidth, rVals[i], 1, forcing, 1);
883 
884  // Note we are using the
885  // LinearAdvectionDiffusionReaction solver here
886  // instead of HelmSolve since m_filtWidth is negative and
887  // so matrices are not positive definite. Ideally
888  // should allow for negative m_filtWidth coefficient in
889  // HelmSolve
890  m_recvField->LinearAdvectionDiffusionReactionSolve(
891  Velocity, forcing, tmpC, -m_filtWidth);
892 
893  m_evalField->BwdTrans(tmpC, field[i]);
894  timer3.Stop();
895 
896  if (m_evalField->GetComm()->GetRank() == 0 &&
897  m_evalField->GetSession()->DefinesCmdLineArgument(
898  "verbose"))
899  {
900  cout << "Smoother time (" << m_recvFieldNames[i]
901  << "): " << timer3.TimePerTest(1) << endl;
902  }
903  }
904  }
905  else
906  {
907  for (int i = 0; i < m_nRecvVars; ++i)
908  {
909  m_recvField->FwdTrans(rVals[i], tmpC);
910  m_evalField->BwdTrans(tmpC, field[i]);
911  }
912  }
913  timer1.Stop();
914 
915  if (m_evalField->GetComm()->GetRank() == 0 &&
916  m_evalField->GetSession()->DefinesCmdLineArgument("verbose"))
917  {
918  cout << "Receive total time: " << timer1.TimePerTest(1) << ", ";
919  cout << "Receive waiting time: " << timer2.TimePerTest(1) << endl;
920  }
921 
922  ReceiveStart();
923  }
924 }
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:225

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

677 {
678  if (m_recvHandle >= 0)
679  {
680  return;
681  }
682 
683  LibUtilities::Timer timer1;
684  timer1.Start();
685  // workaround a bug in cwipi: receiving_field_name should be const char* but
686  // is char*
687  char recFN[10];
688  strcpy(recFN, "dummyName");
689 
690  cwipi_irecv(m_couplingName.c_str(),
691  "ex1",
692  m_recvTag,
693  m_nRecvVars,
694  0,
695  0.0,
696  recFN,
698  &m_recvHandle);
699  timer1.Stop();
700 
701  if (m_evalField->GetComm()->GetRank() == 0 &&
702  m_evalField->GetSession()->DefinesCmdLineArgument("verbose"))
703  {
704  cout << "Receive start time: " << timer1.TimePerTest(1) << endl;
705  }
706 }

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

567 {
568  ASSERTL0(interpField.size() == m_nSendVars, "size mismatch");
569 
570  for (int i = 0; i < m_nSendVars; ++i)
571  {
572  interpField[i] = Array<OneD, NekDouble>(distCoords.size());
573  }
574 
575  if (boost::to_upper_copy(m_config["SENDMETHOD"]) == "NEARESTNEIGHBOUR" ||
576  boost::to_upper_copy(m_config["SENDMETHOD"]) == "SHEPARD")
577  {
578  if (!m_sendInterpolator)
579  {
581  }
582 
585  0, m_sendField);
588  0, interpField);
589  m_sendInterpolator->Interpolate(ptsIn, ptsOut);
590  }
591  else
592  {
593  EvaluateFields(interpField, distCoords);
594  }
595 }
FieldUtils::InterpolatorSharedPtr 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 654 of file CouplingCwipi.cpp.

655 {
656  if (m_sendHandle < 0)
657  {
658  return;
659  }
660 
661  LibUtilities::Timer timer1;
662  timer1.Start();
663  cwipi_wait_issend(m_couplingName.c_str(), m_sendHandle);
664  timer1.Stop();
665 
666  // set to -1 so we dont try finishing a send before a new one was started
667  m_sendHandle = -1;
668 
669  if (m_evalField->GetComm()->GetRank() == 0 &&
670  m_evalField->GetSession()->DefinesCmdLineArgument("verbose"))
671  {
672  cout << "Send waiting time: " << timer1.TimePerTest(1) << endl;
673  }
674 }

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

260 {
261  int oversamp = boost::lexical_cast<int>(m_config["OVERSAMPLE"]);
262 
265  recvGraph->SetExpansionInfoToPointOrder(
266  oversamp + m_evalField->GetExp(0)->GetNumPoints(0));
267 
268  // TODO: DeclareCoeffPhysArrays
269  m_recvField =
271  m_evalField->GetSession(), recvGraph, "DefaultVar");
272 
273  m_oldFields = Array<OneD, Array<OneD, NekDouble> >(m_nRecvVars);
274  m_newFields = Array<OneD, Array<OneD, NekDouble> >(m_nRecvVars);
275  for (int i = 0; i < m_nRecvVars; ++i)
276  {
277  m_oldFields[i] = Array<OneD, NekDouble>(m_evalField->GetTotPoints());
278  m_newFields[i] = Array<OneD, NekDouble>(m_evalField->GetTotPoints());
279  }
280 
281  // define the quadrature points at which we want to receive data
282  m_nPoints = m_recvField->GetTotPoints();
283  Array<OneD, Array<OneD, NekDouble> > coords(3);
284  coords[0] = Array<OneD, NekDouble>(m_nPoints);
285  coords[1] = Array<OneD, NekDouble>(m_nPoints);
286  coords[2] = Array<OneD, NekDouble>(m_nPoints);
287  m_recvField->GetCoords(coords[0], coords[1], coords[2]);
288 
289  m_points = (double *)malloc(sizeof(double) * 3 * m_nPoints);
290  ASSERTL1(m_points != NULL, "malloc failed for m_points");
291 
292  for (int i = 0; i < m_nPoints; ++i)
293  {
294  m_points[3 * i + 0] = double(coords[0][i]);
295 
296  if (m_spacedim > 1)
297  {
298  m_points[3 * i + 1] = double(coords[1][i]);
299  }
300  else
301  {
302  m_points[3 * i + 1] = 0.0;
303  }
304 
305  if (m_spacedim > 2)
306  {
307  m_points[3 * i + 2] = double(coords[2][i]);
308  }
309  else
310  {
311  m_points[3 * i + 2] = 0.0;
312  }
313  }
314 
315  cwipi_set_points_to_locate(m_couplingName.c_str(), m_nPoints, m_points);
316 
317  m_rValsInterl = (double *)malloc(sizeof(double) * m_nPoints * m_nRecvVars);
318  ASSERTL1(m_rValsInterl != NULL, "malloc failed for m_rValsInterl");
319 }
Array< OneD, Array< OneD, NekDouble > > m_oldFields
Array< OneD, Array< OneD, NekDouble > > m_newFields
static MeshGraphSharedPtr Read(const LibUtilities::SessionReaderSharedPtr pSession, LibUtilities::DomainRangeShPtr rng=LibUtilities::NullDomainRangeShPtr, bool fillGraph=true)
Definition: MeshGraph.cpp:113

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

322 {
323  // this array is never used because of our send callback method
324  m_sValsInterl = (double *)malloc(
325  sizeof(double) * m_evalField->GetGraph()->GetNvertices() * m_nSendVars);
326  ASSERTL1(m_sValsInterl != NULL, "malloc failed for m_sValsInterl");
327  for (int i = 0; i < m_evalField->GetGraph()->GetNvertices() * m_nSendVars;
328  ++i)
329  {
330  m_sValsInterl[i] = i;
331  }
332 
333  // register this coupling as sender
334  std::stringstream sst;
335  sst << m_spacedim << "," << m_evalField->GetGraph()->GetNvertices() << ","
336  << m_nSendVars;
337  SendCallbackMap[sst.str()] = std::bind(&CouplingCwipi::SendCallback,
338  this,
339  std::placeholders::_1,
340  std::placeholders::_2);
341  cwipi_set_interpolation_function(m_couplingName.c_str(),
343 }
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 396 of file CouplingCwipi.cpp.

397 {
398  const double *distCoords =
399  cwipi_get_distant_coordinates(m_couplingName.c_str());
400  int nPts = cwipi_get_n_distant_points(m_couplingName.c_str());
401 
402  Array<OneD, Array<OneD, NekDouble> > local(3);
403  for (int i = 0; i < 3; ++i)
404  {
405  local[i] = Array<OneD, NekDouble>(m_evalField->GetTotPoints(), 0.0);
406  }
407  m_evalField->GetCoords(local[0], local[1], local[2]);
410 
411  Array<OneD, Array<OneD, NekDouble> > dist(3);
412  for (int i = 0; i < 3; ++i)
413  {
414  dist[i] = Array<OneD, NekDouble>(nPts);
415  for (int j = 0; j < nPts; ++j)
416  {
417  dist[i][j] = distCoords[3 * j + i];
418  }
419  }
422 
424  if (boost::to_upper_copy(m_config["SENDMETHOD"]) == "SHEPARD")
425  {
426  method = LibUtilities::eShepard;
427  }
430  m_sendInterpolator->CalcWeights(locatPts, distPts);
431  m_sendInterpolator->PrintStatistics();
432 }

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  )
protectedvirtual

Reimplemented from Nektar::SolverUtils::Coupling.

Definition at line 514 of file CouplingCwipi.cpp.

515 {
516  cwipi_delete_coupling(m_couplingName.c_str());
517 }

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

◆ v_Init()

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

Reimplemented from Nektar::SolverUtils::Coupling.

Definition at line 142 of file CouplingCwipi.cpp.

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

Implements Nektar::SolverUtils::Coupling.

Definition at line 708 of file CouplingCwipi.cpp.

712 {
713  if (m_nRecvVars < 1 || m_recvSteps < 1)
714  {
715  return;
716  }
717 
718  Array<OneD, Array<OneD, NekDouble> > recvFields(m_nRecvVars);
719  vector<int> recvVarsToVars =
721  ASSERTL1(m_nRecvVars == recvVarsToVars.size(), "field size mismatch");
722  for (int i = 0; i < recvVarsToVars.size(); ++i)
723  {
724  recvFields[i] = field[recvVarsToVars[i]];
725  }
726 
727  int nq = m_evalField->GetTotPoints();
728 
729  // make sure we have sensible data in old/new recvFields the first time this
730  // method is called
731  if (m_lastReceive < 0)
732  {
733  for (int i = 0; i < m_nRecvVars; ++i)
734  {
735  Vmath::Vcopy(nq, recvFields[i], 1, m_oldFields[i], 1);
736  Vmath::Vcopy(nq, recvFields[i], 1, m_newFields[i], 1);
737  }
738  }
739 
740  if (step >= m_lastReceive + m_recvSteps)
741  {
742  for (int i = 0; i < m_nRecvVars; ++i)
743  {
744  Vmath::Vcopy(nq, m_newFields[i], 1, m_oldFields[i], 1);
745  }
746 
747  ReceiveCwipi(step, time, m_newFields);
748  }
749 
750  NekDouble fact =
752  for (int i = 0; i < m_nRecvVars; ++i)
753  {
754  Vmath::Svtsvtp(nq,
755  fact,
756  m_newFields[i],
757  1,
758  (1 - fact),
759  m_oldFields[i],
760  1,
761  recvFields[i],
762  1);
763  }
764 }
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:135
void Svtsvtp(int n, const T alpha, const T *x, int incx, const T beta, const T *y, int incy, T *z, int incz)
vvtvvtp (scalar times vector plus scalar times vector):
Definition: Vmath.cpp:691
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1199

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 
)
protectedvirtual

Implements Nektar::SolverUtils::Coupling.

Definition at line 597 of file CouplingCwipi.cpp.

602 {
603  if (m_nSendVars < 1 || m_sendSteps < 1)
604  {
605  return;
606  }
607 
608  if (step >= m_lastSend + m_sendSteps)
609  {
610  SendComplete();
611 
612  m_lastSend = step;
613 
614  vector<int> sendVarsToVars =
616  m_sendField = Array<OneD, Array<OneD, NekDouble> >(m_nSendVars);
617  for (int i = 0; i < sendVarsToVars.size(); ++i)
618  {
619  m_sendField[i] = field[sendVarsToVars[i]];
620  }
621 
622  if (m_evalField->GetComm()->GetRank() == 0 &&
623  m_evalField->GetSession()->DefinesCmdLineArgument("verbose"))
624  {
625  cout << "sending fields at i = " << step << ", t = " << time
626  << endl;
627  }
628 
629  LibUtilities::Timer timer1;
630  timer1.Start();
631 
632  char sendFN[10];
633  strcpy(sendFN, "dummyName");
634 
635  cwipi_issend(m_couplingName.c_str(),
636  "ex1",
637  m_sendTag,
638  m_nSendVars,
639  step,
640  time,
641  sendFN,
643  &m_sendHandle);
644  timer1.Stop();
645 
646  if (m_evalField->GetComm()->GetRank() == 0 &&
647  m_evalField->GetSession()->DefinesCmdLineArgument("verbose"))
648  {
649  cout << "Send total time: " << timer1.TimePerTest(1) << endl;
650  }
651  }
652 }
std::vector< std::string > m_sendFieldNames
Definition: Coupling.h:118

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

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

◆ m_connecIdx

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

Definition at line 120 of file CouplingCwipi.h.

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

◆ m_coords

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

Definition at line 119 of file CouplingCwipi.h.

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

◆ m_extrapInterpolator

FieldUtils::InterpolatorSharedPtr Nektar::SolverUtils::CouplingCwipi::m_extrapInterpolator
protected

Definition at line 129 of file CouplingCwipi.h.

Referenced by ExtrapolateFields().

◆ m_filtWidth

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

Definition at line 98 of file CouplingCwipi.h.

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

◆ m_lastReceive

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

Definition at line 114 of file CouplingCwipi.h.

Referenced by ReceiveCwipi(), and v_Receive().

◆ m_lastSend

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

Definition at line 113 of file CouplingCwipi.h.

Referenced by v_Send().

◆ m_newFields

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

Definition at line 105 of file CouplingCwipi.h.

Referenced by SetupReceive(), and v_Receive().

◆ m_nPoints

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

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

Referenced by SetupReceive(), and v_Receive().

◆ m_points

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

Definition at line 118 of file CouplingCwipi.h.

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

◆ m_recvField

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

Definition at line 102 of file CouplingCwipi.h.

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

◆ m_recvHandle

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

Definition at line 108 of file CouplingCwipi.h.

Referenced by ReceiveCwipi(), and ReceiveStart().

◆ m_recvTag

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

Definition at line 111 of file CouplingCwipi.h.

Referenced by ReceiveStart(), and v_Init().

◆ m_rValsInterl

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

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

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

◆ m_sendHandle

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

Definition at line 107 of file CouplingCwipi.h.

Referenced by SendComplete(), and v_Send().

◆ m_sendInterpolator

FieldUtils::InterpolatorSharedPtr Nektar::SolverUtils::CouplingCwipi::m_sendInterpolator
protected

Definition at line 127 of file CouplingCwipi.h.

Referenced by SendCallback(), and SetupSendInterpolation().

◆ m_sendTag

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

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

Referenced by AddElementsToMesh().