Nektar++
Loading...
Searching...
No Matches
CouplingCwipi.cpp
Go to the documentation of this file.
1////////////////////////////////////////////////////////////////////////////////
2//
3// File: CouplingCwipi.cpp
4//
5// For more information, please see: http://www.nektar.info/
6//
7// The MIT License
8//
9// Copyright (c) 2017 Kilian Lackhove
10//
11// Permission is hereby granted, free of charge, to any person obtaining a
12// copy of this software and associated documentation files (the "Software"),
13// to deal in the Software without restriction, including without limitation
14// the rights to use, copy, modify, merge, publish, distribute, sublicense,
15// and/or sell copies of the Software, and to permit persons to whom the
16// Software is furnished to do so, subject to the following conditions:
17//
18// The above copyright notice and this permission notice shall be included
19// in all copies or substantial portions of the Software.
20//
21// THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
22// OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
23// FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
24// THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
25// LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
26// FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
27// DEALINGS IN THE SOFTWARE.
28//
29// Description: CWIPI Exchange class
30//
31////////////////////////////////////////////////////////////////////////////////
32
33#include "CouplingCwipi.h"
34
43
45
46#include <boost/format.hpp>
47
48#define OUTPUT_FREQ 0
49
50namespace Nektar::SolverUtils
51{
52
53std::string CouplingCwipi::className =
54 GetCouplingFactory().RegisterCreatorFunction("Cwipi", CouplingCwipi::create,
55 "Cwipi Coupling");
56
58 const int entities_dim, const int n_local_vertex,
59 [[maybe_unused]] const int n_local_element,
60 [[maybe_unused]] const int n_local_polyhedra, const int n_distant_point,
61 [[maybe_unused]] const double local_coordinates[],
62 [[maybe_unused]] const int local_connectivity_index[],
63 [[maybe_unused]] const int local_connectivity[],
64 [[maybe_unused]] const int local_polyhedra_face_index[],
65 [[maybe_unused]] const int local_polyhedra_cell_to_face_connectivity[],
66 [[maybe_unused]] const int local_polyhedra_face_connectivity_index[],
67 [[maybe_unused]] const int local_polyhedra_face_connectivity[],
68 const double distant_points_coordinates[],
69 [[maybe_unused]] const int distant_points_location[],
70 [[maybe_unused]] const float distant_points_distance[],
71 [[maybe_unused]] const int distant_points_barycentric_coordinates_index[],
72 [[maybe_unused]] const double distant_points_barycentric_coordinates[],
73 const int stride, [[maybe_unused]] const cwipi_solver_type_t solver_type,
74 [[maybe_unused]] const void *local_field, void *distant_field)
75{
76 Array<OneD, Array<OneD, NekDouble>> interpField(stride);
77
78 Array<OneD, Array<OneD, NekDouble>> distCoords(n_distant_point);
79 for (int i = 0; i < n_distant_point; ++i)
80 {
81 distCoords[i] = Array<OneD, NekDouble>(3);
82 for (int j = 0; j < 3; ++j)
83 {
84 distCoords[i][j] = distant_points_coordinates[3 * i + j];
85 }
86 }
87
88 std::stringstream sst;
89 sst << entities_dim << "," << n_local_vertex << "," << stride;
90 SendCallbackMap[sst.str()](interpField, distCoords);
91
92 ASSERTL0(interpField.size() == stride, "size mismatch");
93 ASSERTL0(interpField[0].size() == n_distant_point, "size mismatch");
94
95 for (int i = 0; i < n_distant_point; i++)
96 {
97 for (int j = 0; j < stride; ++j)
98 {
99 ((double *)distant_field)[i * stride + j] = interpField[j][i];
100 }
101 }
102}
103
105 : Coupling(field), m_sendHandle(-1), m_recvHandle(-1), m_lastSend(-1E6),
106 m_lastReceive(-1E6), m_points(nullptr), m_coords(nullptr),
107 m_connecIdx(nullptr), m_connec(nullptr), m_rValsInterl(nullptr),
108 m_sValsInterl(nullptr)
109{
110 // defaults
111 m_config["GEOMTOL"] = "0.1";
112 m_config["LOCALNAME"] = "nektar";
113 m_config["REMOTENAME"] = "precise";
114 m_config["OVERSAMPLE"] = "0";
115 m_config["FILTERWIDTH"] = "-1";
116 m_config["DUMPRAW"] = "0";
117 m_config["SENDMETHOD"] = "EVALUATE";
118 m_config["NOTLOCMETHOD"] = "KEEP";
119}
120
122{
124
125 ReadConfig(m_evalField->GetSession());
126
127 cwipi_add_local_int_control_parameter("nSendVars", m_nSendVars);
128 cwipi_add_local_int_control_parameter("nRecvVars", m_nRecvVars);
129 cwipi_add_local_string_control_parameter(
130 "recvFieldNames", m_config["RECEIVEVARIABLES"].c_str());
131 cwipi_add_local_string_control_parameter("sendFieldNames",
132 m_config["SENDVARIABLES"].c_str());
133
134 // MPI_TAG_UB is guaranteed to be at least 32767, so we make
135 // sure m_recvTag < 32767. Only caveat: m_recvTag is not guaranteed to be
136 // unique.
137 m_recvTag =
138 std::hash<std::string>()(m_couplingName + m_config["REMOTENAME"] +
139 m_config["LOCALNAME"]) %
140 32767;
141
142 cwipi_add_local_int_control_parameter("receiveTag", m_recvTag);
143
144 m_spacedim = m_evalField->GetGraph()->GetSpaceDimension();
145
146 if (m_filtWidth > 0)
147 {
148 m_filtWidth = 2 * M_PI / m_filtWidth;
150 }
151
152 // Init Coupling
153 cwipi_solver_type_t solver_type = CWIPI_SOLVER_CELL_VERTEX;
154 NekDouble geom_tol = std::stod(m_config["GEOMTOL"]);
155 cwipi_create_coupling(
156 m_couplingName.c_str(), CWIPI_COUPLING_PARALLEL_WITH_PARTITIONING,
157 m_config["REMOTENAME"].c_str(), m_spacedim, geom_tol, CWIPI_STATIC_MESH,
158 solver_type, OUTPUT_FREQ, "Ensight Gold", "text");
159 cwipi_synchronize_control_parameter(m_config["REMOTENAME"].c_str());
160
161 if (m_evalField->GetComm()->GetRank() == 0 &&
162 m_evalField->GetSession()->DefinesCmdLineArgument("verbose"))
163 {
164 cwipi_dump_application_properties();
165 }
166
167 m_sendTag = cwipi_get_distant_int_control_parameter(
168 m_config["REMOTENAME"].c_str(), "receiveTag");
169
170 if (cwipi_has_int_parameter(m_config["REMOTENAME"].c_str(), "nRecvVars"))
171 {
172 int remoteNRecvVars = cwipi_get_distant_int_control_parameter(
173 m_config["REMOTENAME"].c_str(), "nRecvVars");
174 ASSERTL0(remoteNRecvVars == m_nSendVars,
175 "Number of local send vars different to remote received vars");
176 }
177
178 if (cwipi_has_int_parameter(m_config["REMOTENAME"].c_str(), "nSendVars"))
179 {
180 int remoteNSendVars = cwipi_get_distant_int_control_parameter(
181 m_config["REMOTENAME"].c_str(), "nSendVars");
182 ASSERTL0(remoteNSendVars == m_nRecvVars,
183 "Number of local receive vars different to remote sent vars");
184 }
185
186 AnnounceMesh();
187
188 if (m_nRecvVars > 0 && m_recvSteps > 0)
189 {
190 SetupReceive();
191 }
192
193 if (m_evalField->GetComm()->GetRank() == 0 &&
194 m_evalField->GetSession()->DefinesCmdLineArgument("verbose"))
195 {
196 std::cout << "locating..." << std::endl;
197 }
198 cwipi_locate(m_couplingName.c_str());
199
200 if (m_nSendVars > 0 && m_sendSteps > 0)
201 {
202 SetupSend();
203 }
204
205 if (m_nRecvVars > 0 && m_recvSteps > 0)
206 {
207 ReceiveStart();
208 }
209}
210
212{
213 free(m_coords);
214 free(m_points);
215 free(m_connec);
216 free(m_connecIdx);
217 free(m_rValsInterl);
218 free(m_sValsInterl);
219}
220
222{
223 ASSERTL0(session->DefinesElement("Nektar/Coupling"),
224 "No Coupling config found");
225
226 m_config["LOCALNAME"] = session->GetCmdLineArgument<std::string>("cwipi");
227
228 TiXmlElement *vCoupling = session->GetElement("Nektar/Coupling");
229 ASSERTL0(vCoupling, "Invalid Coupling config");
230
231 m_filtWidth = std::stod(m_config["FILTERWIDTH"]);
232}
233
235{
236 int oversamp = std::stoi(m_config["OVERSAMPLE"]);
237
240 recvGraph->SetExpansionInfoToPointOrder(
241 oversamp + m_evalField->GetExp(0)->GetNumPoints(0));
242
243 // TODO: DeclareCoeffPhysArrays
245 m_evalField->GetSession(), recvGraph, "DefaultVar");
246
249 for (int i = 0; i < m_nRecvVars; ++i)
250 {
251 m_oldFields[i] = Array<OneD, NekDouble>(m_evalField->GetTotPoints());
252 m_newFields[i] = Array<OneD, NekDouble>(m_evalField->GetTotPoints());
253 }
254
255 // define the quadrature points at which we want to receive data
256 m_nPoints = m_recvField->GetTotPoints();
261 m_recvField->GetCoords(coords[0], coords[1], coords[2]);
262
263 m_points = (double *)malloc(sizeof(double) * 3 * m_nPoints);
264 ASSERTL1(m_points != nullptr, "malloc failed for m_points");
265
266 for (int i = 0; i < m_nPoints; ++i)
267 {
268 m_points[3 * i + 0] = double(coords[0][i]);
269
270 if (m_spacedim > 1)
271 {
272 m_points[3 * i + 1] = double(coords[1][i]);
273 }
274 else
275 {
276 m_points[3 * i + 1] = 0.0;
277 }
278
279 if (m_spacedim > 2)
280 {
281 m_points[3 * i + 2] = double(coords[2][i]);
282 }
283 else
284 {
285 m_points[3 * i + 2] = 0.0;
286 }
287 }
288
289 cwipi_set_points_to_locate(m_couplingName.c_str(), m_nPoints, m_points);
290
291 m_rValsInterl = (double *)malloc(sizeof(double) * m_nPoints * m_nRecvVars);
292 ASSERTL1(m_rValsInterl != nullptr, "malloc failed for m_rValsInterl");
293}
294
296{
297 // this array is never used because of our send callback method
298 m_sValsInterl = (double *)malloc(
299 sizeof(double) * m_evalField->GetGraph()->GetNvertices() * m_nSendVars);
300 ASSERTL1(m_sValsInterl != nullptr, "malloc failed for m_sValsInterl");
301 for (int i = 0; i < m_evalField->GetGraph()->GetNvertices() * m_nSendVars;
302 ++i)
303 {
304 m_sValsInterl[i] = i;
305 }
306
307 // register this coupling as sender
308 std::stringstream sst;
309 sst << m_spacedim << "," << m_evalField->GetGraph()->GetNvertices() << ","
310 << m_nSendVars;
311 SendCallbackMap[sst.str()] =
312 std::bind(&CouplingCwipi::SendCallback, this, std::placeholders::_1,
313 std::placeholders::_2);
314 cwipi_set_interpolation_function(m_couplingName.c_str(),
316}
317
319 Array<OneD, Array<OneD, NekDouble>> interpField,
320 Array<OneD, Array<OneD, NekDouble>> distCoords)
321{
322 int nOutPts = distCoords.size();
323
325 for (int i = 0; i < nOutPts; ++i)
326 {
327 // Obtain Element and LocalCoordinate to interpolate
328 int elmtid = -1;
330 while (elmtid < 0 && tol <= 1E3 * NekConstants::kNekZeroTol)
331 {
332 elmtid = m_evalField->GetExpIndex(distCoords[i], Lcoords, tol);
333 tol *= 2;
334 }
335 if (tol > 2 * NekConstants::kNekZeroTol)
336 {
337 for (int j = 0; j < m_spacedim; ++j)
338 {
339 if (Lcoords[j] < -1 - 0.75 * NekConstants::kNekZeroTol)
340 {
341 Lcoords[j] = -1;
342 }
343 if (Lcoords[j] > 1 + 0.75 * NekConstants::kNekZeroTol)
344 {
345 Lcoords[j] = 1;
346 }
347 }
348 }
349
350 ASSERTL0(elmtid >= 0,
351 "no element found for (" +
352 boost::lexical_cast<std::string>(distCoords[i][0]) + ", " +
353 boost::lexical_cast<std::string>(distCoords[i][1]) + ", " +
354 boost::lexical_cast<std::string>(distCoords[i][2]) + ")");
355
356 int offset = m_evalField->GetPhys_Offset(elmtid);
357
358 for (int f = 0; f < m_nSendVars; ++f)
359 {
360 NekDouble value = m_evalField->GetExp(elmtid)->StdPhysEvaluate(
361 Lcoords, m_sendField[f] + offset);
362
363 ASSERTL0(!(boost::math::isnan)(value), "new value is not a number");
364 interpField[f][i] = value;
365 }
366 }
367}
368
370{
371 const double *distCoords =
372 cwipi_get_distant_coordinates(m_couplingName.c_str());
373 int nPts = cwipi_get_n_distant_points(m_couplingName.c_str());
374
376 for (int i = 0; i < 3; ++i)
377 {
378 local[i] = Array<OneD, NekDouble>(m_evalField->GetTotPoints(), 0.0);
379 }
380 m_evalField->GetCoords(local[0], local[1], local[2]);
383
385 for (int i = 0; i < 3; ++i)
386 {
387 dist[i] = Array<OneD, NekDouble>(nPts);
388 for (int j = 0; j < nPts; ++j)
389 {
390 dist[i][j] = distCoords[3 * j + i];
391 }
392 }
395
397 if (boost::to_upper_copy(m_config["SENDMETHOD"]) == "SHEPARD")
398 {
399 method = LibUtilities::eShepard;
400 }
402 MultiRegions::ExpListSharedPtr>>>::AllocateSharedPtr(method);
403 m_sendInterpolator->CalcWeights(locatPts, distPts);
404 m_sendInterpolator->PrintStatistics();
405}
406
408{
410
411 int nElts = 0, tmp = 0;
412
413 // get Elements
414 if (m_spacedim == 1)
415 {
416 int nSeg = graph->GetGeomMap<SpatialDomains::SegGeom>().size();
417 nElts += nSeg;
418 tmp += nSeg * 2;
419 }
420 else if (m_spacedim == 2)
421 {
422 int nTri = graph->GetGeomMap<SpatialDomains::TriGeom>().size();
423 int nQuad = graph->GetGeomMap<SpatialDomains::QuadGeom>().size();
424 nElts += nTri + nQuad;
425 tmp += nTri * 3 + nQuad * 4;
426 }
427 else if (m_spacedim == 3)
428 {
429 int nTet = graph->GetGeomMap<SpatialDomains::TetGeom>().size();
430 int nPyr = graph->GetGeomMap<SpatialDomains::PyrGeom>().size();
431 int nPrism = graph->GetGeomMap<SpatialDomains::PrismGeom>().size();
432 int nHex = graph->GetGeomMap<SpatialDomains::HexGeom>().size();
433
434 nElts += nTet + nPyr + nPrism + nHex;
435 tmp += 4 * nTet + 5 * nPyr + 6 * nPrism + 8 * nHex;
436 }
437
438 int nVerts = graph->GetNvertices();
439
440 // allocate CWIPI arrays
441 m_coords = (double *)malloc(sizeof(double) * 3 * nVerts);
442 ASSERTL1(m_coords != nullptr, "malloc failed for m_coords");
443 m_connec = (int *)malloc(sizeof(int) * tmp);
444 ASSERTL1(m_connec != nullptr, "malloc failed for m_connec");
445 m_connecIdx = (int *)malloc(sizeof(int) * (nElts + 1));
446 ASSERTL1(m_connecIdx != nullptr, "malloc failed for m_connecIdx");
447
448 m_connecIdx[0] = 0;
449 int coordsPos = 0;
450 int connecPos = 0;
451 int conidxPos = 0;
452
453 if (m_spacedim == 1)
454 {
455 AddElementsToMesh(graph->GetGeomMap<SpatialDomains::SegGeom>(),
456 coordsPos, connecPos, conidxPos);
457 }
458 else if (m_spacedim == 2)
459 {
460 AddElementsToMesh(graph->GetGeomMap<SpatialDomains::TriGeom>(),
461 coordsPos, connecPos, conidxPos);
462 AddElementsToMesh(graph->GetGeomMap<SpatialDomains::QuadGeom>(),
463 coordsPos, connecPos, conidxPos);
464 }
465 else if (m_spacedim == 3)
466 {
467 AddElementsToMesh(graph->GetGeomMap<SpatialDomains::TetGeom>(),
468 coordsPos, connecPos, conidxPos);
469 AddElementsToMesh(graph->GetGeomMap<SpatialDomains::PyrGeom>(),
470 coordsPos, connecPos, conidxPos);
471 AddElementsToMesh(graph->GetGeomMap<SpatialDomains::PrismGeom>(),
472 coordsPos, connecPos, conidxPos);
473 AddElementsToMesh(graph->GetGeomMap<SpatialDomains::HexGeom>(),
474 coordsPos, connecPos, conidxPos);
475 }
476
477 // output the mesh in tecplot format. If this works, CWIPI will be able
478 // to process it, too
479 /*
480 std::cout << "VARIABLES = \"X\", \"Y\", \"Z\", \"U\"" << std::endl;
481 std::cout << "ZONE NODES=" << nVerts << ",ELEMENTS=" << nElts;
482 std::cout << ",DATAPACKING=POINT, ZONETYPE=FETETRAHEDRON" << std::endl;
483 for (int i = 0; i < nVerts; ++i)
484 {
485 std::cout << m_coords[3*i + 0] << " " << m_coords[3*i + 1] << " ";
486 std::cout << m_coords[3*i + 2] << " " << 1.0 << std::endl;
487 }
488 for (int i = 0; i < nElts; ++i)
489 {
490 std::cout << m_connec[i*4 + 0] << " " << m_connec[i*4 + 1] << " ";
491 std::cout << m_connec[i*4 + 2] << " " << m_connec[i*4 + 3] <<
492 std::endl;
493 }
494 */
495
496 cwipi_define_mesh(m_couplingName.c_str(), nVerts, nElts, m_coords,
498}
499
501{
502 cwipi_delete_coupling(m_couplingName.c_str());
503}
504
505template <typename T>
506void CouplingCwipi::AddElementsToMesh(T &geomMap, int &coordsPos,
507 int &connecPos, int &conidxPos)
508{
509 // If map is empty (i.e. no elements), do nothing.
510 if (geomMap.size() == 0)
511 {
512 return;
513 }
514
515 // helper variables
518 int vertID;
519
520 int kNverts = (*geomMap.begin()).second->GetNumVerts();
521
522 // iterate over all elements
523 for (auto [id, geom] : geomMap)
524 {
525 // iterate over the elements vertices
526 for (int j = 0; j < kNverts; ++j)
527 {
528 vert = geom->GetVertex(j);
529 vertID = vert->GetGlobalID();
530
531 // check if we already stored the vertex
532 if (m_vertMap.count(vertID) == 0)
533 {
534 // store the vertex
535 vert->GetCoords(x[0], x[1], x[2]);
536 m_coords[3 * coordsPos + 0] = double(x[0]);
537 m_coords[3 * coordsPos + 1] = double(x[1]);
538 m_coords[3 * coordsPos + 2] = double(x[2]);
539
540 // store the vertex position in the m_coords array
541 m_vertMap[vertID] = coordsPos;
542 coordsPos++;
543 }
544
545 m_connec[connecPos] = m_vertMap[vertID] + 1;
546 connecPos++;
547 }
548
549 m_connecIdx[conidxPos + 1] = m_connecIdx[conidxPos] + kNverts;
550 conidxPos++;
551 }
552}
553
555 Array<OneD, Array<OneD, NekDouble>> &interpField,
556 Array<OneD, Array<OneD, NekDouble>> &distCoords)
557{
558 ASSERTL0(interpField.size() == m_nSendVars, "size mismatch");
559
560 for (int i = 0; i < m_nSendVars; ++i)
561 {
562 interpField[i] = Array<OneD, NekDouble>(distCoords.size());
563 }
564
565 if (boost::to_upper_copy(m_config["SENDMETHOD"]) == "NEARESTNEIGHBOUR" ||
566 boost::to_upper_copy(m_config["SENDMETHOD"]) == "SHEPARD")
567 {
569 {
571 }
572
575 0, m_sendField);
578 0, interpField);
579 m_sendInterpolator->Interpolate(ptsIn, ptsOut);
580 }
581 else
582 {
583 EvaluateFields(interpField, distCoords);
584 }
585}
586
588 const int step, const NekDouble time,
589 const Array<OneD, const Array<OneD, NekDouble>> &field,
590 std::vector<std::string> &varNames)
591{
592 if (m_nSendVars < 1 || m_sendSteps < 1)
593 {
594 return;
595 }
596
597 if (step >= m_lastSend + m_sendSteps)
598 {
599 SendComplete();
600
601 m_lastSend = step;
602
603 std::vector<int> sendVarsToVars =
606 for (int i = 0; i < sendVarsToVars.size(); ++i)
607 {
608 m_sendField[i] = field[sendVarsToVars[i]];
609 }
610
611 if (m_evalField->GetComm()->GetRank() == 0 &&
612 m_evalField->GetSession()->DefinesCmdLineArgument("verbose"))
613 {
614 std::cout << "sending fields at i = " << step << ", t = " << time
615 << std::endl;
616 }
617
618 LibUtilities::Timer timer1;
619 timer1.Start();
620
621 char sendFN[10];
622 strcpy(sendFN, "dummyName");
623
624 cwipi_issend(m_couplingName.c_str(), "ex1", m_sendTag, m_nSendVars,
625 step, time, sendFN, m_sValsInterl, &m_sendHandle);
626 timer1.Stop();
627
628 if (m_evalField->GetComm()->GetRank() == 0 &&
629 m_evalField->GetSession()->DefinesCmdLineArgument("verbose"))
630 {
631 std::cout << "Send total time: " << timer1.TimePerTest(1)
632 << std::endl;
633 }
634 }
635}
636
638{
639 if (m_sendHandle < 0)
640 {
641 return;
642 }
643
644 LibUtilities::Timer timer1;
645 timer1.Start();
646 cwipi_wait_issend(m_couplingName.c_str(), m_sendHandle);
647 timer1.Stop();
648
649 // set to -1 so we dont try finishing a send before a new one was started
650 m_sendHandle = -1;
651
652 if (m_evalField->GetComm()->GetRank() == 0 &&
653 m_evalField->GetSession()->DefinesCmdLineArgument("verbose"))
654 {
655 std::cout << "Send waiting time: " << timer1.TimePerTest(1)
656 << std::endl;
657 }
658}
659
661{
662 if (m_recvHandle >= 0)
663 {
664 return;
665 }
666
667 LibUtilities::Timer timer1;
668 timer1.Start();
669 // workaround a bug in cwipi: receiving_field_name should be const char* but
670 // is char*
671 char recFN[10];
672 strcpy(recFN, "dummyName");
673
674 cwipi_irecv(m_couplingName.c_str(), "ex1", m_recvTag, m_nRecvVars, 0, 0.0,
675 recFN, m_rValsInterl, &m_recvHandle);
676 timer1.Stop();
677
678 if (m_evalField->GetComm()->GetRank() == 0 &&
679 m_evalField->GetSession()->DefinesCmdLineArgument("verbose"))
680 {
681 std::cout << "Receive start time: " << timer1.TimePerTest(1)
682 << std::endl;
683 }
684}
685
686void CouplingCwipi::v_Receive(const int step, const NekDouble time,
688 std::vector<std::string> &varNames)
689{
690 if (m_nRecvVars < 1 || m_recvSteps < 1)
691 {
692 return;
693 }
694
696 std::vector<int> recvVarsToVars =
698 ASSERTL1(m_nRecvVars == recvVarsToVars.size(), "field size mismatch");
699 for (int i = 0; i < recvVarsToVars.size(); ++i)
700 {
701 recvFields[i] = field[recvVarsToVars[i]];
702 }
703
704 int nq = m_evalField->GetTotPoints();
705
706 // make sure we have sensible data in old/new recvFields the first time this
707 // method is called
708 if (m_lastReceive < 0)
709 {
710 for (int i = 0; i < m_nRecvVars; ++i)
711 {
712 Vmath::Vcopy(nq, recvFields[i], 1, m_oldFields[i], 1);
713 Vmath::Vcopy(nq, recvFields[i], 1, m_newFields[i], 1);
714 }
715 }
716
717 if (step >= m_lastReceive + m_recvSteps)
718 {
719 for (int i = 0; i < m_nRecvVars; ++i)
720 {
721 Vmath::Vcopy(nq, m_newFields[i], 1, m_oldFields[i], 1);
722 }
723
724 ReceiveCwipi(step, time, m_newFields);
725 }
726
727 NekDouble fact =
729 for (int i = 0; i < m_nRecvVars; ++i)
730 {
731 Vmath::Svtsvtp(nq, fact, m_newFields[i], 1, (1 - fact), m_oldFields[i],
732 1, recvFields[i], 1);
733 }
734}
735
736void CouplingCwipi::ReceiveCwipi(const int step, const NekDouble time,
738{
739 ASSERTL1(m_nRecvVars == field.size(), "field size mismatch");
740
741 if (m_nRecvVars < 1 || m_recvSteps < 1)
742 {
743 return;
744 }
745
746 if (step >= m_lastReceive + m_recvSteps)
747 {
748 m_lastReceive = step;
749
750 if (m_evalField->GetComm()->GetRank() == 0 &&
751 m_evalField->GetSession()->DefinesCmdLineArgument("verbose"))
752 {
753 std::cout << "waiting for receive at i = " << step
754 << ", t = " << time << std::endl;
755 }
756
757 LibUtilities::Timer timer1, timer2, timer3;
758 timer1.Start();
759
760 Array<OneD, NekDouble> tmpC(m_recvField->GetNcoeffs());
762 for (int i = 0; i < m_nRecvVars; ++i)
763 {
764 rVals[i] = Array<OneD, NekDouble>(m_recvField->GetTotPoints());
765 }
766
767 timer2.Start();
768 cwipi_wait_irecv(m_couplingName.c_str(), m_recvHandle);
769 timer2.Stop();
770
771 // set to -1 so we know we can start receiving again
772 m_recvHandle = -1;
773
774 int nNotLoc = cwipi_get_n_not_located_points(m_couplingName.c_str());
775 Array<OneD, int> notLoc;
776 if (nNotLoc != 0)
777 {
778 std::cout << "WARNING: relocating " << nNotLoc << " of "
779 << m_nPoints << " points" << std::endl;
780
781 const int *tmp =
782 cwipi_get_not_located_points(m_couplingName.c_str());
783 notLoc = Array<OneD, int>(nNotLoc);
784 for (int i = 0; i < nNotLoc; ++i)
785 {
786 notLoc[i] = tmp[i] - 1;
787 }
788
789 if (boost::to_upper_copy(m_config["NOTLOCMETHOD"]) == "KEEP")
790 {
791 // interpolate from m_evalField to m_recvField
792 for (int i = 0; i < m_nRecvVars; ++i)
793 {
794 m_evalField->FwdTrans(field[i], tmpC);
795 m_recvField->BwdTrans(tmpC, rVals[i]);
796 }
797 }
798 }
799
800 for (int i = 0, locPos = 0, intPos = 0; i < m_nPoints; ++i)
801 {
802 if (locPos < nNotLoc && notLoc[locPos] == i)
803 {
804 // keep the original value of field[j][i]
805 locPos++;
806 }
807 else
808 {
809 for (int j = 0; j < m_nRecvVars; ++j)
810 {
811 rVals[j][i] = m_rValsInterl[intPos * m_nRecvVars + j];
812 }
813 intPos++;
814 }
815 }
816
817 if (boost::to_upper_copy(m_config["NOTLOCMETHOD"]) == "EXTRAPOLATE")
818 {
819 int doExtrapolate = 0;
820 if (nNotLoc != 0)
821 {
822 doExtrapolate = 1;
823 }
824 m_evalField->GetSession()->GetComm()->AllReduce(
825 doExtrapolate, LibUtilities::ReduceMax);
826 if (doExtrapolate > 0)
827 {
828 ExtrapolateFields(rVals, notLoc);
829 }
830 }
831
832 if (m_config["DUMPRAW"] != "0")
833 {
834 DumpRawFields(time, rVals);
835 }
836
837 if (m_filtWidth > 0)
838 {
839 for (int i = 0; i < m_nRecvVars; ++i)
840 {
841 timer3.Start();
842
844
846 for (int j = 0; j < m_spacedim; ++j)
847 {
848 Velocity[j] = Array<OneD, NekDouble>(m_nPoints, 0.0);
849 }
850
851 Vmath::Smul(m_nPoints, -m_filtWidth, rVals[i], 1, forcing, 1);
852
854 StdRegions::VarCoeffMap varcoeffs;
855 MultiRegions::VarFactorsMap varfactors =
857
859
860 // Set advection velocities
861 StdRegions::VarCoeffType varcoefftypes[] = {
864 for (int j = 0; j < m_spacedim; j++)
865 {
866 varcoeffs[varcoefftypes[j]] = Velocity[j];
867 }
868
869 // Note we are using the
870 // LinearAdvectionDiffusionReaction solver here
871 // instead of HelmSolve since m_filtWidth is negative and
872 // so matrices are not positive definite. Ideally
873 // should allow for negative m_filtWidth coefficient in
874 // HelmSolve
875 // m_recvField->LinearAdvectionDiffusionReactionSolve(
876 // Velocity, forcing, tmpC, -m_filtWidth);
877 m_recvField->LinearAdvectionDiffusionReactionSolve(
878 forcing, tmpC, factors, varcoeffs, varfactors);
879
880 m_evalField->BwdTrans(tmpC, field[i]);
881 timer3.Stop();
882
883 if (m_evalField->GetComm()->GetRank() == 0 &&
884 m_evalField->GetSession()->DefinesCmdLineArgument(
885 "verbose"))
886 {
887 std::cout << "Smoother time (" << m_recvFieldNames[i]
888 << "): " << timer3.TimePerTest(1) << std::endl;
889 }
890 }
891 }
892 else
893 {
894 for (int i = 0; i < m_nRecvVars; ++i)
895 {
896 m_recvField->FwdTrans(rVals[i], tmpC);
897 m_evalField->BwdTrans(tmpC, field[i]);
898 }
899 }
900 timer1.Stop();
901
902 if (m_evalField->GetComm()->GetRank() == 0 &&
903 m_evalField->GetSession()->DefinesCmdLineArgument("verbose"))
904 {
905 std::cout << "Receive total time: " << timer1.TimePerTest(1)
906 << ", ";
907 std::cout << "Receive waiting time: " << timer2.TimePerTest(1)
908 << std::endl;
909 }
910
911 ReceiveStart();
912 }
913}
914
917{
918 LibUtilities::Timer timer1, timer2;
919 timer1.Start();
920
921 int totvars = 3 + m_nRecvVars;
922 int nNotLoc = notLoc.size();
923 int nranks = m_evalField->GetSession()->GetComm()->GetSize();
924
925 Array<OneD, int> thisNLoc(1, m_nPoints - nNotLoc);
926 Array<OneD, int> sizeMap(nranks);
927 m_evalField->GetSession()->GetComm()->AllGather(thisNLoc, sizeMap);
928
929 Array<OneD, int> offsetMap(nranks);
930 offsetMap[0] = 0;
931 for (int i = 1; i < nranks; ++i)
932 {
933 offsetMap[i] = offsetMap[i - 1] + sizeMap[i - 1];
934 }
935 int totNLoc = offsetMap[nranks - 1] + sizeMap[nranks - 1];
936
937 Array<OneD, Array<OneD, NekDouble>> allVals(totvars);
938 for (int i = 0; i < 3; ++i)
939 {
940 allVals[i] = Array<OneD, NekDouble>(m_nPoints);
941 }
942 m_recvField->GetCoords(allVals[0], allVals[1], allVals[2]);
943
944 if (m_spacedim < 3)
945 {
946 Vmath::Zero(m_nPoints, allVals[2], 1);
947 }
948 if (m_spacedim < 2)
949 {
950 Vmath::Zero(m_nPoints, allVals[1], 1);
951 }
952
953 for (int i = 0; i < m_nRecvVars; ++i)
954 {
955 allVals[3 + i] = rVals[i];
956 }
957
958 Array<OneD, Array<OneD, NekDouble>> locatedVals(totvars);
959 for (int i = 0; i < totvars; ++i)
960 {
961 locatedVals[i] = Array<OneD, NekDouble>(totNLoc, -42.0);
962 }
963
964 // only copy points from allVals to locatedVals that were located
965 int offset = offsetMap[m_evalField->GetSession()->GetComm()->GetRank()];
966 for (int i = 0, intPos = 0, locPos = 0; i < m_nPoints; ++i)
967 {
968 if (locPos < nNotLoc && notLoc[locPos] == i)
969 {
970 // do nothing
971 locPos++;
972 }
973 else
974 {
975 for (int k = 0; k < totvars; ++k)
976 {
977 locatedVals[k][offset + intPos] = allVals[k][i];
978 }
979 intPos++;
980 }
981 }
982
983 // send all located points to all ranks. This is probably horribly expensive
984 timer2.Start();
985 for (int i = 0; i < totvars; ++i)
986 {
987 m_evalField->GetSession()->GetComm()->AllGatherv(locatedVals[i],
988 sizeMap, offsetMap);
989 }
990 timer2.Stop();
991
992 if (nNotLoc > 0)
993 {
996 3, locatedVals);
997
998 Array<OneD, Array<OneD, NekDouble>> notLocVals(totvars);
999 for (int j = 0; j < totvars; ++j)
1000 {
1001 notLocVals[j] = Array<OneD, NekDouble>(nNotLoc);
1002 for (int i = 0; i < nNotLoc; ++i)
1003 {
1004 notLocVals[j][i] = allVals[j][notLoc[i]];
1005 }
1006 }
1009 3, notLocVals);
1010
1011 // perform a nearest neighbour interpolation from locatedVals to the not
1012 // located rVals
1014 {
1016 std::vector<MultiRegions::ExpListSharedPtr>>>::
1017 AllocateSharedPtr(LibUtilities::eNearestNeighbour);
1018 m_extrapInterpolator->CalcWeights(locatedPts, notlocPts);
1019 m_extrapInterpolator->PrintStatistics();
1020 }
1021 m_extrapInterpolator->Interpolate(locatedPts, notlocPts);
1022
1023 for (int j = 3; j < totvars; ++j)
1024 {
1025 for (int i = 0; i < nNotLoc; ++i)
1026 {
1027 allVals[j][notLoc[i]] = notLocVals[j][i];
1028 }
1029 }
1030 }
1031
1032 timer1.Stop();
1033 if (m_evalField->GetComm()->GetRank() == 0 &&
1034 m_evalField->GetSession()->DefinesCmdLineArgument("verbose"))
1035 {
1036 std::cout << "ExtrapolateFields total time: " << timer1.TimePerTest(1);
1037 std::cout << " (AllGathervI: " << timer2.TimePerTest(1) << ")"
1038 << std::endl;
1039 }
1040}
1041
1044{
1045 LibUtilities::Timer timer1;
1046 timer1.Start();
1047
1048#if (defined _WIN32 && _MSC_VER < 1900)
1049 // We need this to make sure boost::format has always
1050 // two digits in the exponents of Scientific notation.
1051 unsigned int old_exponent_format;
1052 old_exponent_format = _set_output_format(_TWO_DIGIT_EXPONENT);
1053 filename = boost::str(boost::format(filename) % m_time);
1054 _set_output_format(old_exponent_format);
1055#else
1056 std::string filename =
1057 boost::str(boost::format(m_config["DUMPRAW"]) % time);
1058#endif
1059
1061 for (int i = 0; i < 3; ++i)
1062 {
1063 tmp[i] = Array<OneD, NekDouble>(m_recvField->GetTotPoints(), 0.0);
1064 }
1065 m_recvField->GetCoords(tmp[0], tmp[1], tmp[2]);
1066
1067 for (int i = 0; i < m_nRecvVars; ++i)
1068 {
1069 tmp[m_spacedim + i] = rVals[i];
1070 }
1071
1072 LibUtilities::CsvIO csvIO(m_evalField->GetSession()->GetComm());
1076 csvIO.Write(filename, rvPts);
1077
1078 timer1.Stop();
1079 if (m_evalField->GetComm()->GetRank() == 0 &&
1080 m_evalField->GetSession()->DefinesCmdLineArgument("verbose"))
1081 {
1082 std::cout << "DumpRawFields total time: " << timer1.TimePerTest(1)
1083 << std::endl;
1084 }
1085}
1086} // namespace Nektar::SolverUtils
#define OUTPUT_FREQ
#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.
Definition CsvIO.cpp:68
NekDouble TimePerTest(unsigned int n)
Returns amount of seconds per iteration in a test with n iterations.
Definition Timer.cpp:65
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)
Array< OneD, Array< OneD, NekDouble > > m_oldFields
SOLVER_UTILS_EXPORT void v_Receive(const int step, const NekDouble time, Array< OneD, Array< OneD, NekDouble > > &field, std::vector< std::string > &varNames) override
Array< OneD, Array< OneD, NekDouble > > m_sendField
SOLVER_UTILS_EXPORT CouplingCwipi(MultiRegions::ExpListSharedPtr field)
SOLVER_UTILS_EXPORT void v_Init() override
void ExtrapolateFields(Array< OneD, Array< OneD, NekDouble > > &rVals, Array< OneD, int > &notLoc)
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)
SOLVER_UTILS_EXPORT void v_Finalize() override
SOLVER_UTILS_EXPORT ~CouplingCwipi() override
static CouplingSharedPtr create(MultiRegions::ExpListSharedPtr field)
Creates an instance of this class.
Array< OneD, Array< OneD, NekDouble > > m_newFields
void AddElementsToMesh(T &geomMap, int &coordsPos, int &connecPos, int &conidxPos)
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
void DumpRawFields(const NekDouble time, Array< OneD, Array< OneD, NekDouble > > rVals)
MultiRegions::ExpListSharedPtr m_recvField
virtual SOLVER_UTILS_EXPORT void v_Init()
Definition Coupling.cpp:57
MultiRegions::ExpListSharedPtr m_evalField
Definition Coupling.h:107
std::vector< std::string > m_sendFieldNames
Definition Coupling.h:110
CouplingConfigMap m_config
Definition Coupling.h:105
std::vector< std::string > m_recvFieldNames
Definition Coupling.h:114
SOLVER_UTILS_EXPORT std::vector< int > GenerateVariableMapping(std::vector< std::string > &vars, std::vector< std::string > &transVars)
Definition Coupling.cpp:131
int GetGlobalID(void) const
Get the ID of this object.
Definition Geometry.h:322
PointGeom * GetVertex(int i) const
Returns vertex i of this object.
Definition Geometry.h:361
int GetNumVerts() const
Get the number of vertices of this object.
Definition Geometry.h:403
static MeshGraphSharedPtr Read(const LibUtilities::SessionReaderSharedPtr pSession, LibUtilities::DomainRangeShPtr rng=LibUtilities::NullDomainRangeShPtr, bool fillGraph=true, SpatialDomains::MeshGraphSharedPtr partitionedGraph=nullptr)
void GetCoords(NekDouble &x, NekDouble &y, NekDouble &z)
Definition PointGeom.cpp:69
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< PtsField > PtsFieldSharedPtr
Definition PtsField.h:184
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.
Definition Coupling.cpp:40
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition MeshGraph.h:217
std::map< ConstFactorType, NekDouble > ConstFactorMap
std::map< StdRegions::VarCoeffType, VarCoeffEntry > VarCoeffMap
void Svtsvtp(int n, const T alpha, const T *x, int incx, const T beta, const T *y, int incy, T *z, int incz)
Svtsvtp (scalar times vector plus scalar times vector):
Definition Vmath.hpp:473
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition Vmath.hpp:100
void Zero(int n, T *x, const int incx)
Zero vector.
Definition Vmath.hpp:273
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition Vmath.hpp:825