45#include <boost/core/ignore_unused.hpp>
46#include <boost/format.hpp>
47#include <boost/functional/hash.hpp>
63 const int entities_dim,
const int n_local_vertex,
const int n_local_element,
64 const int n_local_polyhedra,
const int n_distant_point,
65 const double local_coordinates[],
const int local_connectivity_index[],
66 const int local_connectivity[],
const int local_polyhedra_face_index[],
67 const int local_polyhedra_cell_to_face_connectivity[],
68 const int local_polyhedra_face_connectivity_index[],
69 const int local_polyhedra_face_connectivity[],
70 const double distant_points_coordinates[],
71 const int distant_points_location[],
const float distant_points_distance[],
72 const int distant_points_barycentric_coordinates_index[],
73 const double distant_points_barycentric_coordinates[],
const int stride,
74 const cwipi_solver_type_t solver_type,
const void *local_field,
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);
89 for (
int i = 0; i < n_distant_point; ++i)
92 for (
int j = 0; j < 3; ++j)
94 distCoords[i][j] = distant_points_coordinates[3 * i + j];
98 std::stringstream sst;
99 sst << entities_dim <<
"," << n_local_vertex <<
"," << stride;
102 ASSERTL0(interpField.size() == stride,
"size mismatch");
103 ASSERTL0(interpField[0].size() == n_distant_point,
"size mismatch");
105 for (
int i = 0; i < n_distant_point; i++)
107 for (
int j = 0; j < stride; ++j)
109 ((
double *)distant_field)[i * stride + j] = interpField[j][i];
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)
126 m_config[
"SENDMETHOD"] =
"EVALUATE";
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",
151 cwipi_add_local_int_control_parameter(
"receiveTag",
m_recvTag);
162 cwipi_solver_type_t solver_type = CWIPI_SOLVER_CELL_VERTEX;
164 cwipi_create_coupling(
165 m_couplingName.c_str(), CWIPI_COUPLING_PARALLEL_WITH_PARTITIONING,
168 cwipi_synchronize_control_parameter(
m_config[
"REMOTENAME"].c_str());
171 m_evalField->GetSession()->DefinesCmdLineArgument(
"verbose"))
173 cwipi_dump_application_properties();
176 m_sendTag = cwipi_get_distant_int_control_parameter(
177 m_config[
"REMOTENAME"].c_str(),
"receiveTag");
179 if (cwipi_has_int_parameter(
m_config[
"REMOTENAME"].c_str(),
"nRecvVars"))
181 int remoteNRecvVars = cwipi_get_distant_int_control_parameter(
182 m_config[
"REMOTENAME"].c_str(),
"nRecvVars");
184 "Number of local send vars different to remote received vars");
187 if (cwipi_has_int_parameter(
m_config[
"REMOTENAME"].c_str(),
"nSendVars"))
189 int remoteNSendVars = cwipi_get_distant_int_control_parameter(
190 m_config[
"REMOTENAME"].c_str(),
"nSendVars");
192 "Number of local receive vars different to remote sent vars");
203 m_evalField->GetSession()->DefinesCmdLineArgument(
"verbose"))
205 cout <<
"locating..." << endl;
232 ASSERTL0(session->DefinesElement(
"Nektar/Coupling"),
233 "No Coupling config found");
235 m_config[
"LOCALNAME"] = session->GetCmdLineArgument<std::string>(
"cwipi");
237 TiXmlElement *vCoupling = session->GetElement(
"Nektar/Coupling");
238 ASSERTL0(vCoupling,
"Invalid Coupling config");
245 int oversamp = boost::lexical_cast<int>(
m_config[
"OVERSAMPLE"]);
249 recvGraph->SetExpansionInfoToPointOrder(
250 oversamp +
m_evalField->GetExp(0)->GetNumPoints(0));
254 m_evalField->GetSession(), recvGraph,
"DefaultVar");
270 m_recvField->GetCoords(coords[0], coords[1], coords[2]);
277 m_points[3 * i + 0] = double(coords[0][i]);
281 m_points[3 * i + 1] = double(coords[1][i]);
290 m_points[3 * i + 2] = double(coords[2][i]);
317 std::stringstream sst;
322 std::placeholders::_2);
331 int nOutPts = distCoords.size();
334 for (
int i = 0; i < nOutPts; ++i)
341 elmtid =
m_evalField->GetExpIndex(distCoords[i], Lcoords, tol);
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]) +
")");
372 ASSERTL0(!(boost::math::isnan)(value),
"new value is not a number");
373 interpField[f][i] = value;
380 const double *distCoords =
385 for (
int i = 0; i < 3; ++i)
389 m_evalField->GetCoords(local[0], local[1], local[2]);
394 for (
int i = 0; i < 3; ++i)
397 for (
int j = 0; j < nPts; ++j)
399 dist[i][j] = distCoords[3 * j + i];
406 if (boost::to_upper_copy(
m_config[
"SENDMETHOD"]) ==
"SHEPARD")
430 seggeom = graph->GetAllSegGeoms();
434 trigeom = graph->GetAllTriGeoms();
435 quadgeom = graph->GetAllQuadGeoms();
439 tetgeom = graph->GetAllTetGeoms();
440 pyrgeom = graph->GetAllPyrGeoms();
441 prismgeom = graph->GetAllPrismGeoms();
442 hexgeom = graph->GetAllHexGeoms();
445 int nVerts = graph->GetNvertices();
446 int nElts = seggeom.size() + trigeom.size() + quadgeom.size() +
447 tetgeom.size() + pyrgeom.size() + prismgeom.size() +
451 m_coords = (
double *)malloc(
sizeof(
double) * 3 * nVerts);
453 int tmp = 2 * seggeom.size() + 3 * trigeom.size() + 4 * quadgeom.size() +
454 4 * tetgeom.size() + 5 * pyrgeom.size() + 6 * prismgeom.size() +
456 m_connec = (
int *)malloc(
sizeof(
int) * tmp);
458 m_connecIdx = (
int *)malloc(
sizeof(
int) * (nElts + 1));
510 int kNverts = T::mapped_type::element_type::kNverts;
513 for (
auto it = geom.begin(); it != geom.end(); it++)
516 for (
int j = 0; j < kNverts; ++j)
518 vert = it->second->GetVertex(j);
519 vertID = vert->GetVid();
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]);
555 if (boost::to_upper_copy(
m_config[
"SENDMETHOD"]) ==
"NEARESTNEIGHBOUR" ||
556 boost::to_upper_copy(
m_config[
"SENDMETHOD"]) ==
"SHEPARD")
580 vector<string> &varNames)
593 vector<int> sendVarsToVars =
596 for (
int i = 0; i < sendVarsToVars.size(); ++i)
602 m_evalField->GetSession()->DefinesCmdLineArgument(
"verbose"))
604 cout <<
"sending fields at i = " << step <<
", t = " << time
612 strcpy(sendFN,
"dummyName");
619 m_evalField->GetSession()->DefinesCmdLineArgument(
"verbose"))
621 cout <<
"Send total time: " << timer1.
TimePerTest(1) << endl;
642 m_evalField->GetSession()->DefinesCmdLineArgument(
"verbose"))
644 cout <<
"Send waiting time: " << timer1.
TimePerTest(1) << endl;
660 strcpy(recFN,
"dummyName");
667 m_evalField->GetSession()->DefinesCmdLineArgument(
"verbose"))
669 cout <<
"Receive start time: " << timer1.
TimePerTest(1) << endl;
675 vector<string> &varNames)
683 vector<int> recvVarsToVars =
686 for (
int i = 0; i < recvVarsToVars.size(); ++i)
688 recvFields[i] = field[recvVarsToVars[i]];
719 1, recvFields[i], 1);
738 m_evalField->GetSession()->DefinesCmdLineArgument(
"verbose"))
740 cout <<
"waiting for receive at i = " << step <<
", t = " << time
761 int nNotLoc = cwipi_get_n_not_located_points(
m_couplingName.c_str());
765 cout <<
"WARNING: relocating " << nNotLoc <<
" of " <<
m_nPoints
766 <<
" points" << endl;
771 for (
int i = 0; i < nNotLoc; ++i)
773 notLoc[i] = tmp[i] - 1;
776 if (boost::to_upper_copy(
m_config[
"NOTLOCMETHOD"]) ==
"KEEP")
787 for (
int i = 0, locPos = 0, intPos = 0; i <
m_nPoints; ++i)
789 if (locPos < nNotLoc && notLoc[locPos] == i)
804 if (boost::to_upper_copy(
m_config[
"NOTLOCMETHOD"]) ==
"EXTRAPOLATE")
806 int doExtrapolate = 0;
813 if (doExtrapolate > 0)
853 varcoeffs[varcoefftypes[j]] = Velocity[j];
864 m_recvField->LinearAdvectionDiffusionReactionSolve(
865 forcing, tmpC,
factors, varcoeffs, varfactors);
890 m_evalField->GetSession()->DefinesCmdLineArgument(
"verbose"))
892 cout <<
"Receive total time: " << timer1.
TimePerTest(1) <<
", ";
893 cout <<
"Receive waiting time: " << timer2.
TimePerTest(1) << endl;
907 int nNotLoc = notLoc.size();
908 int nranks =
m_evalField->GetSession()->GetComm()->GetSize();
912 m_evalField->GetSession()->GetComm()->AllGather(thisNLoc, sizeMap);
916 for (
int i = 1; i < nranks; ++i)
918 offsetMap[i] = offsetMap[i - 1] + sizeMap[i - 1];
920 int totNLoc = offsetMap[nranks - 1] + sizeMap[nranks - 1];
923 for (
int i = 0; i < 3; ++i)
927 m_recvField->GetCoords(allVals[0], allVals[1], allVals[2]);
940 allVals[3 + i] = rVals[i];
944 for (
int i = 0; i < totvars; ++i)
950 int offset = offsetMap[
m_evalField->GetSession()->GetComm()->GetRank()];
951 for (
int i = 0, intPos = 0, locPos = 0; i <
m_nPoints; ++i)
953 if (locPos < nNotLoc && notLoc[locPos] == i)
960 for (
int k = 0; k < totvars; ++k)
962 locatedVals[k][offset + intPos] = allVals[k][i];
970 for (
int i = 0; i < totvars; ++i)
972 m_evalField->GetSession()->GetComm()->AllGatherv(locatedVals[i],
984 for (
int j = 0; j < totvars; ++j)
987 for (
int i = 0; i < nNotLoc; ++i)
989 notLocVals[j][i] = allVals[j][notLoc[i]];
1001 std::vector<MultiRegions::ExpListSharedPtr>>>::
1008 for (
int j = 3; j < totvars; ++j)
1010 for (
int i = 0; i < nNotLoc; ++i)
1012 allVals[j][notLoc[i]] = notLocVals[j][i];
1019 m_evalField->GetSession()->DefinesCmdLineArgument(
"verbose"))
1021 cout <<
"ExtrapolateFields total time: " << timer1.
TimePerTest(1);
1022 cout <<
" (AllGathervI: " << timer2.
TimePerTest(1) <<
")" << endl;
1032#if (defined _WIN32 && _MSC_VER < 1900)
1035 unsigned int old_exponent_format;
1036 old_exponent_format = _set_output_format(_TWO_DIGIT_EXPONENT);
1038 _set_output_format(old_exponent_format);
1040 std::string filename =
1045 for (
int i = 0; i < 3; ++i)
1060 csvIO.
Write(filename, rvPts);
1064 m_evalField->GetSession()->DefinesCmdLineArgument(
"verbose"))
1066 cout <<
"DumpRawFields total time: " << timer1.
TimePerTest(1) << endl;
#define ASSERTL0(condition, msg)
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
A class that contains algorithms for interpolation between pts fields, expansions and different meshe...
void Write(const std::string &outFile, const PtsFieldSharedPtr &ptsField, const bool backup=false)
Save a pts field to a file.
NekDouble TimePerTest(unsigned int n)
Returns amount of seconds per iteration in a test with n iterations.
General purpose memory allocation routines with the ability to allocate from thread specific memory p...
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
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 std::string className
Array< OneD, Array< OneD, NekDouble > > m_oldFields
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
void AddElementsToMesh(T geom, int &coordsPos, int &connecPos, int &conidxPos)
Array< OneD, Array< OneD, NekDouble > > m_sendField
virtual SOLVER_UTILS_EXPORT ~CouplingCwipi()
SOLVER_UTILS_EXPORT CouplingCwipi(MultiRegions::ExpListSharedPtr field)
virtual SOLVER_UTILS_EXPORT void v_Init() override
void ExtrapolateFields(Array< OneD, Array< OneD, NekDouble > > &rVals, Array< OneD, int > ¬Loc)
std::shared_ptr< FieldUtils::Interpolator< std::vector< MultiRegions::ExpListSharedPtr > > > m_sendInterpolator
void EvaluateFields(Array< OneD, Array< OneD, NekDouble > > interpField, Array< OneD, Array< OneD, NekDouble > > distCoords)
virtual SOLVER_UTILS_EXPORT void v_Finalize() override
static CouplingSharedPtr create(MultiRegions::ExpListSharedPtr field)
Creates an instance of this class.
Array< OneD, Array< OneD, NekDouble > > m_newFields
void SetupSendInterpolation()
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
SOLVER_UTILS_EXPORT void SendCallback(Array< OneD, Array< OneD, NekDouble > > &interpField, Array< OneD, Array< OneD, NekDouble > > &distCoords)
void ReadConfig(LibUtilities::SessionReaderSharedPtr session)
void ReceiveCwipi(const int step, const NekDouble time, Array< OneD, Array< OneD, NekDouble > > &field)
std::shared_ptr< FieldUtils::Interpolator< std::vector< MultiRegions::ExpListSharedPtr > > > m_extrapInterpolator
std::map< int, int > m_vertMap
void DumpRawFields(const NekDouble time, Array< OneD, Array< OneD, NekDouble > > rVals)
MultiRegions::ExpListSharedPtr m_recvField
virtual SOLVER_UTILS_EXPORT void v_Init()
MultiRegions::ExpListSharedPtr m_evalField
std::string m_couplingName
std::vector< std::string > m_sendFieldNames
CouplingConfigMap m_config
std::vector< std::string > m_recvFieldNames
SOLVER_UTILS_EXPORT std::vector< int > GenerateVariableMapping(std::vector< std::string > &vars, std::vector< std::string > &transVars)
static MeshGraphSharedPtr Read(const LibUtilities::SessionReaderSharedPtr pSession, LibUtilities::DomainRangeShPtr rng=LibUtilities::NullDomainRangeShPtr, bool fillGraph=true, SpatialDomains::MeshGraphSharedPtr partitionedGraph=nullptr)
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< PtsField > PtsFieldSharedPtr
static VarFactorsMap NullVarFactorsMap
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
std::map< StdRegions::ConstFactorType, Array< OneD, NekDouble > > VarFactorsMap
static const NekDouble kNekZeroTol
static std::map< std::string, SendCallbackType > SendCallbackMap
CouplingFactory & GetCouplingFactory()
Declaration of the Coupling factory singleton.
std::map< int, TriGeomSharedPtr > TriGeomMap
std::map< int, PyrGeomSharedPtr > PyrGeomMap
std::map< int, QuadGeomSharedPtr > QuadGeomMap
std::map< int, SegGeomSharedPtr > SegGeomMap
std::map< int, TetGeomSharedPtr > TetGeomMap
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
std::map< int, PrismGeomSharedPtr > PrismGeomMap
std::shared_ptr< PointGeom > PointGeomSharedPtr
std::map< int, HexGeomSharedPtr > HexGeomMap
std::map< ConstFactorType, NekDouble > ConstFactorMap
std::map< StdRegions::VarCoeffType, VarCoeffEntry > VarCoeffMap
StdRegions::ConstFactorMap factors
The above copyright notice and this permission notice shall be included.
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):
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
void Zero(int n, T *x, const int incx)
Zero vector.
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)