Nektar++
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#include <boost/functional/hash.hpp>
48
49#define OUTPUT_FREQ 0
50
51namespace Nektar::SolverUtils
52{
53
54using namespace std;
55
56std::string CouplingCwipi::className =
57 GetCouplingFactory().RegisterCreatorFunction("Cwipi", CouplingCwipi::create,
58 "Cwipi Coupling");
59
61 const int entities_dim, const int n_local_vertex,
62 [[maybe_unused]] const int n_local_element,
63 [[maybe_unused]] const int n_local_polyhedra, const int n_distant_point,
64 [[maybe_unused]] const double local_coordinates[],
65 [[maybe_unused]] const int local_connectivity_index[],
66 [[maybe_unused]] const int local_connectivity[],
67 [[maybe_unused]] const int local_polyhedra_face_index[],
68 [[maybe_unused]] const int local_polyhedra_cell_to_face_connectivity[],
69 [[maybe_unused]] const int local_polyhedra_face_connectivity_index[],
70 [[maybe_unused]] const int local_polyhedra_face_connectivity[],
71 const double distant_points_coordinates[],
72 [[maybe_unused]] const int distant_points_location[],
73 [[maybe_unused]] const float distant_points_distance[],
74 [[maybe_unused]] const int distant_points_barycentric_coordinates_index[],
75 [[maybe_unused]] const double distant_points_barycentric_coordinates[],
76 const int stride, [[maybe_unused]] const cwipi_solver_type_t solver_type,
77 [[maybe_unused]] const void *local_field, void *distant_field)
78{
79 Array<OneD, Array<OneD, NekDouble>> interpField(stride);
80
81 Array<OneD, Array<OneD, NekDouble>> distCoords(n_distant_point);
82 for (int i = 0; i < n_distant_point; ++i)
83 {
84 distCoords[i] = Array<OneD, NekDouble>(3);
85 for (int j = 0; j < 3; ++j)
86 {
87 distCoords[i][j] = distant_points_coordinates[3 * i + j];
88 }
89 }
90
91 std::stringstream sst;
92 sst << entities_dim << "," << n_local_vertex << "," << stride;
93 SendCallbackMap[sst.str()](interpField, distCoords);
94
95 ASSERTL0(interpField.size() == stride, "size mismatch");
96 ASSERTL0(interpField[0].size() == n_distant_point, "size mismatch");
97
98 for (int i = 0; i < n_distant_point; i++)
99 {
100 for (int j = 0; j < stride; ++j)
101 {
102 ((double *)distant_field)[i * stride + j] = interpField[j][i];
103 }
104 }
105}
106
108 : Coupling(field), m_sendHandle(-1), m_recvHandle(-1), m_lastSend(-1E6),
109 m_lastReceive(-1E6), m_points(nullptr), m_coords(nullptr),
110 m_connecIdx(nullptr), m_connec(nullptr), m_rValsInterl(nullptr),
111 m_sValsInterl(nullptr)
112{
113 // defaults
114 m_config["GEOMTOL"] = "0.1";
115 m_config["LOCALNAME"] = "nektar";
116 m_config["REMOTENAME"] = "precise";
117 m_config["OVERSAMPLE"] = "0";
118 m_config["FILTERWIDTH"] = "-1";
119 m_config["DUMPRAW"] = "0";
120 m_config["SENDMETHOD"] = "EVALUATE";
121 m_config["NOTLOCMETHOD"] = "KEEP";
122}
123
125{
127
128 ReadConfig(m_evalField->GetSession());
129
130 cwipi_add_local_int_control_parameter("nSendVars", m_nSendVars);
131 cwipi_add_local_int_control_parameter("nRecvVars", m_nRecvVars);
132 cwipi_add_local_string_control_parameter(
133 "recvFieldNames", m_config["RECEIVEVARIABLES"].c_str());
134 cwipi_add_local_string_control_parameter("sendFieldNames",
135 m_config["SENDVARIABLES"].c_str());
136
137 // MPI_TAG_UB is guaranteed to be at least 32767, so we make
138 // sure m_recvTag < 32767. Only caveat: m_recvTag is not guaranteed to be
139 // unique.
140 m_recvTag =
141 boost::hash<std::string>()(m_couplingName + m_config["REMOTENAME"] +
142 m_config["LOCALNAME"]) %
143 32767;
144
145 cwipi_add_local_int_control_parameter("receiveTag", m_recvTag);
146
147 m_spacedim = m_evalField->GetGraph()->GetSpaceDimension();
148
149 if (m_filtWidth > 0)
150 {
151 m_filtWidth = 2 * M_PI / m_filtWidth;
153 }
154
155 // Init Coupling
156 cwipi_solver_type_t solver_type = CWIPI_SOLVER_CELL_VERTEX;
157 NekDouble geom_tol = std::stod(m_config["GEOMTOL"]);
158 cwipi_create_coupling(
159 m_couplingName.c_str(), CWIPI_COUPLING_PARALLEL_WITH_PARTITIONING,
160 m_config["REMOTENAME"].c_str(), m_spacedim, geom_tol, CWIPI_STATIC_MESH,
161 solver_type, OUTPUT_FREQ, "Ensight Gold", "text");
162 cwipi_synchronize_control_parameter(m_config["REMOTENAME"].c_str());
163
164 if (m_evalField->GetComm()->GetRank() == 0 &&
165 m_evalField->GetSession()->DefinesCmdLineArgument("verbose"))
166 {
167 cwipi_dump_application_properties();
168 }
169
170 m_sendTag = cwipi_get_distant_int_control_parameter(
171 m_config["REMOTENAME"].c_str(), "receiveTag");
172
173 if (cwipi_has_int_parameter(m_config["REMOTENAME"].c_str(), "nRecvVars"))
174 {
175 int remoteNRecvVars = cwipi_get_distant_int_control_parameter(
176 m_config["REMOTENAME"].c_str(), "nRecvVars");
177 ASSERTL0(remoteNRecvVars == m_nSendVars,
178 "Number of local send vars different to remote received vars");
179 }
180
181 if (cwipi_has_int_parameter(m_config["REMOTENAME"].c_str(), "nSendVars"))
182 {
183 int remoteNSendVars = cwipi_get_distant_int_control_parameter(
184 m_config["REMOTENAME"].c_str(), "nSendVars");
185 ASSERTL0(remoteNSendVars == m_nRecvVars,
186 "Number of local receive vars different to remote sent vars");
187 }
188
189 AnnounceMesh();
190
191 if (m_nRecvVars > 0 && m_recvSteps > 0)
192 {
193 SetupReceive();
194 }
195
196 if (m_evalField->GetComm()->GetRank() == 0 &&
197 m_evalField->GetSession()->DefinesCmdLineArgument("verbose"))
198 {
199 cout << "locating..." << endl;
200 }
201 cwipi_locate(m_couplingName.c_str());
202
203 if (m_nSendVars > 0 && m_sendSteps > 0)
204 {
205 SetupSend();
206 }
207
208 if (m_nRecvVars > 0 && m_recvSteps > 0)
209 {
210 ReceiveStart();
211 }
212}
213
215{
216 free(m_coords);
217 free(m_points);
218 free(m_connec);
219 free(m_connecIdx);
220 free(m_rValsInterl);
221 free(m_sValsInterl);
222}
223
225{
226 ASSERTL0(session->DefinesElement("Nektar/Coupling"),
227 "No Coupling config found");
228
229 m_config["LOCALNAME"] = session->GetCmdLineArgument<std::string>("cwipi");
230
231 TiXmlElement *vCoupling = session->GetElement("Nektar/Coupling");
232 ASSERTL0(vCoupling, "Invalid Coupling config");
233
234 m_filtWidth = std::stod(m_config["FILTERWIDTH"]);
235}
236
238{
239 int oversamp = boost::lexical_cast<int>(m_config["OVERSAMPLE"]);
240
243 recvGraph->SetExpansionInfoToPointOrder(
244 oversamp + m_evalField->GetExp(0)->GetNumPoints(0));
245
246 // TODO: DeclareCoeffPhysArrays
248 m_evalField->GetSession(), recvGraph, "DefaultVar");
249
252 for (int i = 0; i < m_nRecvVars; ++i)
253 {
254 m_oldFields[i] = Array<OneD, NekDouble>(m_evalField->GetTotPoints());
255 m_newFields[i] = Array<OneD, NekDouble>(m_evalField->GetTotPoints());
256 }
257
258 // define the quadrature points at which we want to receive data
259 m_nPoints = m_recvField->GetTotPoints();
264 m_recvField->GetCoords(coords[0], coords[1], coords[2]);
265
266 m_points = (double *)malloc(sizeof(double) * 3 * m_nPoints);
267 ASSERTL1(m_points != nullptr, "malloc failed for m_points");
268
269 for (int i = 0; i < m_nPoints; ++i)
270 {
271 m_points[3 * i + 0] = double(coords[0][i]);
272
273 if (m_spacedim > 1)
274 {
275 m_points[3 * i + 1] = double(coords[1][i]);
276 }
277 else
278 {
279 m_points[3 * i + 1] = 0.0;
280 }
281
282 if (m_spacedim > 2)
283 {
284 m_points[3 * i + 2] = double(coords[2][i]);
285 }
286 else
287 {
288 m_points[3 * i + 2] = 0.0;
289 }
290 }
291
292 cwipi_set_points_to_locate(m_couplingName.c_str(), m_nPoints, m_points);
293
294 m_rValsInterl = (double *)malloc(sizeof(double) * m_nPoints * m_nRecvVars);
295 ASSERTL1(m_rValsInterl != nullptr, "malloc failed for m_rValsInterl");
296}
297
299{
300 // this array is never used because of our send callback method
301 m_sValsInterl = (double *)malloc(
302 sizeof(double) * m_evalField->GetGraph()->GetNvertices() * m_nSendVars);
303 ASSERTL1(m_sValsInterl != nullptr, "malloc failed for m_sValsInterl");
304 for (int i = 0; i < m_evalField->GetGraph()->GetNvertices() * m_nSendVars;
305 ++i)
306 {
307 m_sValsInterl[i] = i;
308 }
309
310 // register this coupling as sender
311 std::stringstream sst;
312 sst << m_spacedim << "," << m_evalField->GetGraph()->GetNvertices() << ","
313 << m_nSendVars;
314 SendCallbackMap[sst.str()] =
315 std::bind(&CouplingCwipi::SendCallback, this, std::placeholders::_1,
316 std::placeholders::_2);
317 cwipi_set_interpolation_function(m_couplingName.c_str(),
319}
320
322 Array<OneD, Array<OneD, NekDouble>> interpField,
323 Array<OneD, Array<OneD, NekDouble>> distCoords)
324{
325 int nOutPts = distCoords.size();
326
328 for (int i = 0; i < nOutPts; ++i)
329 {
330 // Obtain Element and LocalCoordinate to interpolate
331 int elmtid = -1;
333 while (elmtid < 0 && tol <= 1E3 * NekConstants::kNekZeroTol)
334 {
335 elmtid = m_evalField->GetExpIndex(distCoords[i], Lcoords, tol);
336 tol *= 2;
337 }
338 if (tol > 2 * NekConstants::kNekZeroTol)
339 {
340 for (int j = 0; j < m_spacedim; ++j)
341 {
342 if (Lcoords[j] < -1 - 0.75 * NekConstants::kNekZeroTol)
343 {
344 Lcoords[j] = -1;
345 }
346 if (Lcoords[j] > 1 + 0.75 * NekConstants::kNekZeroTol)
347 {
348 Lcoords[j] = 1;
349 }
350 }
351 }
352
353 ASSERTL0(elmtid >= 0,
354 "no element found for (" +
355 boost::lexical_cast<string>(distCoords[i][0]) + ", " +
356 boost::lexical_cast<string>(distCoords[i][1]) + ", " +
357 boost::lexical_cast<string>(distCoords[i][2]) + ")");
358
359 int offset = m_evalField->GetPhys_Offset(elmtid);
360
361 for (int f = 0; f < m_nSendVars; ++f)
362 {
363 NekDouble value = m_evalField->GetExp(elmtid)->StdPhysEvaluate(
364 Lcoords, m_sendField[f] + offset);
365
366 ASSERTL0(!(boost::math::isnan)(value), "new value is not a number");
367 interpField[f][i] = value;
368 }
369 }
370}
371
373{
374 const double *distCoords =
375 cwipi_get_distant_coordinates(m_couplingName.c_str());
376 int nPts = cwipi_get_n_distant_points(m_couplingName.c_str());
377
379 for (int i = 0; i < 3; ++i)
380 {
381 local[i] = Array<OneD, NekDouble>(m_evalField->GetTotPoints(), 0.0);
382 }
383 m_evalField->GetCoords(local[0], local[1], local[2]);
386
388 for (int i = 0; i < 3; ++i)
389 {
390 dist[i] = Array<OneD, NekDouble>(nPts);
391 for (int j = 0; j < nPts; ++j)
392 {
393 dist[i][j] = distCoords[3 * j + i];
394 }
395 }
398
400 if (boost::to_upper_copy(m_config["SENDMETHOD"]) == "SHEPARD")
401 {
402 method = LibUtilities::eShepard;
403 }
405 MultiRegions::ExpListSharedPtr>>>::AllocateSharedPtr(method);
406 m_sendInterpolator->CalcWeights(locatPts, distPts);
407 m_sendInterpolator->PrintStatistics();
408}
409
411{
413
414 // get Elements
422 if (m_spacedim == 1)
423 {
424 seggeom = graph->GetAllSegGeoms();
425 }
426 else if (m_spacedim == 2)
427 {
428 trigeom = graph->GetAllTriGeoms();
429 quadgeom = graph->GetAllQuadGeoms();
430 }
431 else if (m_spacedim == 3)
432 {
433 tetgeom = graph->GetAllTetGeoms();
434 pyrgeom = graph->GetAllPyrGeoms();
435 prismgeom = graph->GetAllPrismGeoms();
436 hexgeom = graph->GetAllHexGeoms();
437 };
438
439 int nVerts = graph->GetNvertices();
440 int nElts = seggeom.size() + trigeom.size() + quadgeom.size() +
441 tetgeom.size() + pyrgeom.size() + prismgeom.size() +
442 hexgeom.size();
443
444 // allocate CWIPI arrays
445 m_coords = (double *)malloc(sizeof(double) * 3 * nVerts);
446 ASSERTL1(m_coords != nullptr, "malloc failed for m_coords");
447 int tmp = 2 * seggeom.size() + 3 * trigeom.size() + 4 * quadgeom.size() +
448 4 * tetgeom.size() + 5 * pyrgeom.size() + 6 * prismgeom.size() +
449 8 * hexgeom.size();
450 m_connec = (int *)malloc(sizeof(int) * tmp);
451 ASSERTL1(m_connec != nullptr, "malloc failed for m_connec");
452 m_connecIdx = (int *)malloc(sizeof(int) * (nElts + 1));
453 ASSERTL1(m_connecIdx != nullptr, "malloc failed for m_connecIdx");
454
455 m_connecIdx[0] = 0;
456 int coordsPos = 0;
457 int connecPos = 0;
458 int conidxPos = 0;
459
460 AddElementsToMesh(seggeom, coordsPos, connecPos, conidxPos);
461 AddElementsToMesh(trigeom, coordsPos, connecPos, conidxPos);
462 AddElementsToMesh(quadgeom, coordsPos, connecPos, conidxPos);
463 AddElementsToMesh(tetgeom, coordsPos, connecPos, conidxPos);
464 AddElementsToMesh(pyrgeom, coordsPos, connecPos, conidxPos);
465 AddElementsToMesh(prismgeom, coordsPos, connecPos, conidxPos);
466 AddElementsToMesh(hexgeom, coordsPos, connecPos, conidxPos);
467
468 // output the mesh in tecplot format. If this works, CWIPI will be able
469 // to process it, too
470 /*
471 cout << "VARIABLES = \"X\", \"Y\", \"Z\", \"U\"" << endl;
472 cout << "ZONE NODES=" << nVerts << ",ELEMENTS=" << nElts;
473 cout << ",DATAPACKING=POINT, ZONETYPE=FETETRAHEDRON" << endl;
474 for (int i = 0; i < nVerts; ++i)
475 {
476 cout << m_coords[3*i + 0] << " " << m_coords[3*i + 1] << " ";
477 cout << m_coords[3*i + 2] << " " << 1.0 << endl;
478 }
479 for (int i = 0; i < nElts; ++i)
480 {
481 cout << m_connec[i*4 + 0] << " " << m_connec[i*4 + 1] << " ";
482 cout << m_connec[i*4 + 2] << " " << m_connec[i*4 + 3] << endl;
483 }
484 */
485
486 cwipi_define_mesh(m_couplingName.c_str(), nVerts, nElts, m_coords,
488}
489
491{
492 cwipi_delete_coupling(m_couplingName.c_str());
493}
494
495template <typename T>
496void CouplingCwipi::AddElementsToMesh(T geom, int &coordsPos, int &connecPos,
497 int &conidxPos)
498{
499 // helper variables
502 int vertID;
503
504 int kNverts = T::mapped_type::element_type::kNverts;
505
506 // iterate over all elements
507 for (auto it = geom.begin(); it != geom.end(); it++)
508 {
509 // iterate over the elements vertices
510 for (int j = 0; j < kNverts; ++j)
511 {
512 vert = it->second->GetVertex(j);
513 vertID = vert->GetVid();
514
515 // check if we already stored the vertex
516 if (m_vertMap.count(vertID) == 0)
517 {
518 // store the vertex
519 vert->GetCoords(x[0], x[1], x[2]);
520 m_coords[3 * coordsPos + 0] = double(x[0]);
521 m_coords[3 * coordsPos + 1] = double(x[1]);
522 m_coords[3 * coordsPos + 2] = double(x[2]);
523
524 // store the vertex position in the m_coords array
525 m_vertMap[vertID] = coordsPos;
526 coordsPos++;
527 }
528
529 m_connec[connecPos] = m_vertMap[vertID] + 1;
530 connecPos++;
531 }
532
533 m_connecIdx[conidxPos + 1] = m_connecIdx[conidxPos] + kNverts;
534 conidxPos++;
535 }
536}
537
539 Array<OneD, Array<OneD, NekDouble>> &interpField,
540 Array<OneD, Array<OneD, NekDouble>> &distCoords)
541{
542 ASSERTL0(interpField.size() == m_nSendVars, "size mismatch");
543
544 for (int i = 0; i < m_nSendVars; ++i)
545 {
546 interpField[i] = Array<OneD, NekDouble>(distCoords.size());
547 }
548
549 if (boost::to_upper_copy(m_config["SENDMETHOD"]) == "NEARESTNEIGHBOUR" ||
550 boost::to_upper_copy(m_config["SENDMETHOD"]) == "SHEPARD")
551 {
553 {
555 }
556
559 0, m_sendField);
562 0, interpField);
563 m_sendInterpolator->Interpolate(ptsIn, ptsOut);
564 }
565 else
566 {
567 EvaluateFields(interpField, distCoords);
568 }
569}
570
572 const int step, const NekDouble time,
573 const Array<OneD, const Array<OneD, NekDouble>> &field,
574 vector<string> &varNames)
575{
576 if (m_nSendVars < 1 || m_sendSteps < 1)
577 {
578 return;
579 }
580
581 if (step >= m_lastSend + m_sendSteps)
582 {
583 SendComplete();
584
585 m_lastSend = step;
586
587 vector<int> sendVarsToVars =
590 for (int i = 0; i < sendVarsToVars.size(); ++i)
591 {
592 m_sendField[i] = field[sendVarsToVars[i]];
593 }
594
595 if (m_evalField->GetComm()->GetRank() == 0 &&
596 m_evalField->GetSession()->DefinesCmdLineArgument("verbose"))
597 {
598 cout << "sending fields at i = " << step << ", t = " << time
599 << endl;
600 }
601
602 LibUtilities::Timer timer1;
603 timer1.Start();
604
605 char sendFN[10];
606 strcpy(sendFN, "dummyName");
607
608 cwipi_issend(m_couplingName.c_str(), "ex1", m_sendTag, m_nSendVars,
609 step, time, sendFN, m_sValsInterl, &m_sendHandle);
610 timer1.Stop();
611
612 if (m_evalField->GetComm()->GetRank() == 0 &&
613 m_evalField->GetSession()->DefinesCmdLineArgument("verbose"))
614 {
615 cout << "Send total time: " << timer1.TimePerTest(1) << endl;
616 }
617 }
618}
619
621{
622 if (m_sendHandle < 0)
623 {
624 return;
625 }
626
627 LibUtilities::Timer timer1;
628 timer1.Start();
629 cwipi_wait_issend(m_couplingName.c_str(), m_sendHandle);
630 timer1.Stop();
631
632 // set to -1 so we dont try finishing a send before a new one was started
633 m_sendHandle = -1;
634
635 if (m_evalField->GetComm()->GetRank() == 0 &&
636 m_evalField->GetSession()->DefinesCmdLineArgument("verbose"))
637 {
638 cout << "Send waiting time: " << timer1.TimePerTest(1) << endl;
639 }
640}
641
643{
644 if (m_recvHandle >= 0)
645 {
646 return;
647 }
648
649 LibUtilities::Timer timer1;
650 timer1.Start();
651 // workaround a bug in cwipi: receiving_field_name should be const char* but
652 // is char*
653 char recFN[10];
654 strcpy(recFN, "dummyName");
655
656 cwipi_irecv(m_couplingName.c_str(), "ex1", m_recvTag, m_nRecvVars, 0, 0.0,
657 recFN, m_rValsInterl, &m_recvHandle);
658 timer1.Stop();
659
660 if (m_evalField->GetComm()->GetRank() == 0 &&
661 m_evalField->GetSession()->DefinesCmdLineArgument("verbose"))
662 {
663 cout << "Receive start time: " << timer1.TimePerTest(1) << endl;
664 }
665}
666
667void CouplingCwipi::v_Receive(const int step, const NekDouble time,
669 vector<string> &varNames)
670{
671 if (m_nRecvVars < 1 || m_recvSteps < 1)
672 {
673 return;
674 }
675
677 vector<int> recvVarsToVars =
679 ASSERTL1(m_nRecvVars == recvVarsToVars.size(), "field size mismatch");
680 for (int i = 0; i < recvVarsToVars.size(); ++i)
681 {
682 recvFields[i] = field[recvVarsToVars[i]];
683 }
684
685 int nq = m_evalField->GetTotPoints();
686
687 // make sure we have sensible data in old/new recvFields the first time this
688 // method is called
689 if (m_lastReceive < 0)
690 {
691 for (int i = 0; i < m_nRecvVars; ++i)
692 {
693 Vmath::Vcopy(nq, recvFields[i], 1, m_oldFields[i], 1);
694 Vmath::Vcopy(nq, recvFields[i], 1, m_newFields[i], 1);
695 }
696 }
697
698 if (step >= m_lastReceive + m_recvSteps)
699 {
700 for (int i = 0; i < m_nRecvVars; ++i)
701 {
702 Vmath::Vcopy(nq, m_newFields[i], 1, m_oldFields[i], 1);
703 }
704
705 ReceiveCwipi(step, time, m_newFields);
706 }
707
708 NekDouble fact =
710 for (int i = 0; i < m_nRecvVars; ++i)
711 {
712 Vmath::Svtsvtp(nq, fact, m_newFields[i], 1, (1 - fact), m_oldFields[i],
713 1, recvFields[i], 1);
714 }
715}
716
717void CouplingCwipi::ReceiveCwipi(const int step, const NekDouble time,
719{
720 ASSERTL1(m_nRecvVars == field.size(), "field size mismatch");
721
722 if (m_nRecvVars < 1 || m_recvSteps < 1)
723 {
724 return;
725 }
726
727 if (step >= m_lastReceive + m_recvSteps)
728 {
729 m_lastReceive = step;
730
731 if (m_evalField->GetComm()->GetRank() == 0 &&
732 m_evalField->GetSession()->DefinesCmdLineArgument("verbose"))
733 {
734 cout << "waiting for receive at i = " << step << ", t = " << time
735 << endl;
736 }
737
738 LibUtilities::Timer timer1, timer2, timer3;
739 timer1.Start();
740
741 Array<OneD, NekDouble> tmpC(m_recvField->GetNcoeffs());
743 for (int i = 0; i < m_nRecvVars; ++i)
744 {
745 rVals[i] = Array<OneD, NekDouble>(m_recvField->GetTotPoints());
746 }
747
748 timer2.Start();
749 cwipi_wait_irecv(m_couplingName.c_str(), m_recvHandle);
750 timer2.Stop();
751
752 // set to -1 so we know we can start receiving again
753 m_recvHandle = -1;
754
755 int nNotLoc = cwipi_get_n_not_located_points(m_couplingName.c_str());
756 Array<OneD, int> notLoc;
757 if (nNotLoc != 0)
758 {
759 cout << "WARNING: relocating " << nNotLoc << " of " << m_nPoints
760 << " points" << endl;
761
762 const int *tmp =
763 cwipi_get_not_located_points(m_couplingName.c_str());
764 notLoc = Array<OneD, int>(nNotLoc);
765 for (int i = 0; i < nNotLoc; ++i)
766 {
767 notLoc[i] = tmp[i] - 1;
768 }
769
770 if (boost::to_upper_copy(m_config["NOTLOCMETHOD"]) == "KEEP")
771 {
772 // interpolate from m_evalField to m_recvField
773 for (int i = 0; i < m_nRecvVars; ++i)
774 {
775 m_evalField->FwdTrans(field[i], tmpC);
776 m_recvField->BwdTrans(tmpC, rVals[i]);
777 }
778 }
779 }
780
781 for (int i = 0, locPos = 0, intPos = 0; i < m_nPoints; ++i)
782 {
783 if (locPos < nNotLoc && notLoc[locPos] == i)
784 {
785 // keep the original value of field[j][i]
786 locPos++;
787 }
788 else
789 {
790 for (int j = 0; j < m_nRecvVars; ++j)
791 {
792 rVals[j][i] = m_rValsInterl[intPos * m_nRecvVars + j];
793 }
794 intPos++;
795 }
796 }
797
798 if (boost::to_upper_copy(m_config["NOTLOCMETHOD"]) == "EXTRAPOLATE")
799 {
800 int doExtrapolate = 0;
801 if (nNotLoc != 0)
802 {
803 doExtrapolate = 1;
804 }
805 m_evalField->GetSession()->GetComm()->AllReduce(
806 doExtrapolate, LibUtilities::ReduceMax);
807 if (doExtrapolate > 0)
808 {
809 ExtrapolateFields(rVals, notLoc);
810 }
811 }
812
813 if (m_config["DUMPRAW"] != "0")
814 {
815 DumpRawFields(time, rVals);
816 }
817
818 if (m_filtWidth > 0)
819 {
820 for (int i = 0; i < m_nRecvVars; ++i)
821 {
822 timer3.Start();
823
825
827 for (int j = 0; j < m_spacedim; ++j)
828 {
829 Velocity[j] = Array<OneD, NekDouble>(m_nPoints, 0.0);
830 }
831
832 Vmath::Smul(m_nPoints, -m_filtWidth, rVals[i], 1, forcing, 1);
833
835 StdRegions::VarCoeffMap varcoeffs;
836 MultiRegions::VarFactorsMap varfactors =
838
840
841 // Set advection velocities
842 StdRegions::VarCoeffType varcoefftypes[] = {
845 for (int j = 0; j < m_spacedim; j++)
846 {
847 varcoeffs[varcoefftypes[j]] = Velocity[j];
848 }
849
850 // Note we are using the
851 // LinearAdvectionDiffusionReaction solver here
852 // instead of HelmSolve since m_filtWidth is negative and
853 // so matrices are not positive definite. Ideally
854 // should allow for negative m_filtWidth coefficient in
855 // HelmSolve
856 // m_recvField->LinearAdvectionDiffusionReactionSolve(
857 // Velocity, forcing, tmpC, -m_filtWidth);
858 m_recvField->LinearAdvectionDiffusionReactionSolve(
859 forcing, tmpC, factors, varcoeffs, varfactors);
860
861 m_evalField->BwdTrans(tmpC, field[i]);
862 timer3.Stop();
863
864 if (m_evalField->GetComm()->GetRank() == 0 &&
865 m_evalField->GetSession()->DefinesCmdLineArgument(
866 "verbose"))
867 {
868 cout << "Smoother time (" << m_recvFieldNames[i]
869 << "): " << timer3.TimePerTest(1) << endl;
870 }
871 }
872 }
873 else
874 {
875 for (int i = 0; i < m_nRecvVars; ++i)
876 {
877 m_recvField->FwdTrans(rVals[i], tmpC);
878 m_evalField->BwdTrans(tmpC, field[i]);
879 }
880 }
881 timer1.Stop();
882
883 if (m_evalField->GetComm()->GetRank() == 0 &&
884 m_evalField->GetSession()->DefinesCmdLineArgument("verbose"))
885 {
886 cout << "Receive total time: " << timer1.TimePerTest(1) << ", ";
887 cout << "Receive waiting time: " << timer2.TimePerTest(1) << endl;
888 }
889
890 ReceiveStart();
891 }
892}
893
896{
897 LibUtilities::Timer timer1, timer2;
898 timer1.Start();
899
900 int totvars = 3 + m_nRecvVars;
901 int nNotLoc = notLoc.size();
902 int nranks = m_evalField->GetSession()->GetComm()->GetSize();
903
904 Array<OneD, int> thisNLoc(1, m_nPoints - nNotLoc);
905 Array<OneD, int> sizeMap(nranks);
906 m_evalField->GetSession()->GetComm()->AllGather(thisNLoc, sizeMap);
907
908 Array<OneD, int> offsetMap(nranks);
909 offsetMap[0] = 0;
910 for (int i = 1; i < nranks; ++i)
911 {
912 offsetMap[i] = offsetMap[i - 1] + sizeMap[i - 1];
913 }
914 int totNLoc = offsetMap[nranks - 1] + sizeMap[nranks - 1];
915
916 Array<OneD, Array<OneD, NekDouble>> allVals(totvars);
917 for (int i = 0; i < 3; ++i)
918 {
919 allVals[i] = Array<OneD, NekDouble>(m_nPoints);
920 }
921 m_recvField->GetCoords(allVals[0], allVals[1], allVals[2]);
922
923 if (m_spacedim < 3)
924 {
925 Vmath::Zero(m_nPoints, allVals[2], 1);
926 }
927 if (m_spacedim < 2)
928 {
929 Vmath::Zero(m_nPoints, allVals[1], 1);
930 }
931
932 for (int i = 0; i < m_nRecvVars; ++i)
933 {
934 allVals[3 + i] = rVals[i];
935 }
936
937 Array<OneD, Array<OneD, NekDouble>> locatedVals(totvars);
938 for (int i = 0; i < totvars; ++i)
939 {
940 locatedVals[i] = Array<OneD, NekDouble>(totNLoc, -42.0);
941 }
942
943 // only copy points from allVals to locatedVals that were located
944 int offset = offsetMap[m_evalField->GetSession()->GetComm()->GetRank()];
945 for (int i = 0, intPos = 0, locPos = 0; i < m_nPoints; ++i)
946 {
947 if (locPos < nNotLoc && notLoc[locPos] == i)
948 {
949 // do nothing
950 locPos++;
951 }
952 else
953 {
954 for (int k = 0; k < totvars; ++k)
955 {
956 locatedVals[k][offset + intPos] = allVals[k][i];
957 }
958 intPos++;
959 }
960 }
961
962 // send all located points to all ranks. This is probably horribly expensive
963 timer2.Start();
964 for (int i = 0; i < totvars; ++i)
965 {
966 m_evalField->GetSession()->GetComm()->AllGatherv(locatedVals[i],
967 sizeMap, offsetMap);
968 }
969 timer2.Stop();
970
971 if (nNotLoc > 0)
972 {
975 3, locatedVals);
976
977 Array<OneD, Array<OneD, NekDouble>> notLocVals(totvars);
978 for (int j = 0; j < totvars; ++j)
979 {
980 notLocVals[j] = Array<OneD, NekDouble>(nNotLoc);
981 for (int i = 0; i < nNotLoc; ++i)
982 {
983 notLocVals[j][i] = allVals[j][notLoc[i]];
984 }
985 }
988 3, notLocVals);
989
990 // perform a nearest neighbour interpolation from locatedVals to the not
991 // located rVals
993 {
995 std::vector<MultiRegions::ExpListSharedPtr>>>::
996 AllocateSharedPtr(LibUtilities::eNearestNeighbour);
997 m_extrapInterpolator->CalcWeights(locatedPts, notlocPts);
998 m_extrapInterpolator->PrintStatistics();
999 }
1000 m_extrapInterpolator->Interpolate(locatedPts, notlocPts);
1001
1002 for (int j = 3; j < totvars; ++j)
1003 {
1004 for (int i = 0; i < nNotLoc; ++i)
1005 {
1006 allVals[j][notLoc[i]] = notLocVals[j][i];
1007 }
1008 }
1009 }
1010
1011 timer1.Stop();
1012 if (m_evalField->GetComm()->GetRank() == 0 &&
1013 m_evalField->GetSession()->DefinesCmdLineArgument("verbose"))
1014 {
1015 cout << "ExtrapolateFields total time: " << timer1.TimePerTest(1);
1016 cout << " (AllGathervI: " << timer2.TimePerTest(1) << ")" << endl;
1017 }
1018}
1019
1022{
1023 LibUtilities::Timer timer1;
1024 timer1.Start();
1025
1026#if (defined _WIN32 && _MSC_VER < 1900)
1027 // We need this to make sure boost::format has always
1028 // two digits in the exponents of Scientific notation.
1029 unsigned int old_exponent_format;
1030 old_exponent_format = _set_output_format(_TWO_DIGIT_EXPONENT);
1031 filename = boost::str(boost::format(filename) % m_time);
1032 _set_output_format(old_exponent_format);
1033#else
1034 std::string filename =
1035 boost::str(boost::format(m_config["DUMPRAW"]) % time);
1036#endif
1037
1039 for (int i = 0; i < 3; ++i)
1040 {
1041 tmp[i] = Array<OneD, NekDouble>(m_recvField->GetTotPoints(), 0.0);
1042 }
1043 m_recvField->GetCoords(tmp[0], tmp[1], tmp[2]);
1044
1045 for (int i = 0; i < m_nRecvVars; ++i)
1046 {
1047 tmp[m_spacedim + i] = rVals[i];
1048 }
1049
1050 LibUtilities::CsvIO csvIO(m_evalField->GetSession()->GetComm());
1054 csvIO.Write(filename, rvPts);
1055
1056 timer1.Stop();
1057 if (m_evalField->GetComm()->GetRank() == 0 &&
1058 m_evalField->GetSession()->DefinesCmdLineArgument("verbose"))
1059 {
1060 cout << "DumpRawFields total time: " << timer1.TimePerTest(1) << endl;
1061 }
1062}
1063} // namespace Nektar::SolverUtils
#define OUTPUT_FREQ
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:242
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
Definition: CouplingCwipi.h:96
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
Definition: CouplingCwipi.h:92
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.
Definition: CouplingCwipi.h:56
Array< OneD, Array< OneD, NekDouble > > m_newFields
Definition: CouplingCwipi.h:97
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
Definition: CouplingCwipi.h:94
virtual SOLVER_UTILS_EXPORT void v_Init()
Definition: Coupling.cpp:59
MultiRegions::ExpListSharedPtr m_evalField
Definition: Coupling.h:111
std::vector< std::string > m_sendFieldNames
Definition: Coupling.h:114
CouplingConfigMap m_config
Definition: Coupling.h:109
std::vector< std::string > m_recvFieldNames
Definition: Coupling.h:118
SOLVER_UTILS_EXPORT std::vector< int > GenerateVariableMapping(std::vector< std::string > &vars, std::vector< std::string > &transVars)
Definition: Coupling.cpp:132
static MeshGraphSharedPtr Read(const LibUtilities::SessionReaderSharedPtr pSession, LibUtilities::DomainRangeShPtr rng=LibUtilities::NullDomainRangeShPtr, bool fillGraph=true, SpatialDomains::MeshGraphSharedPtr partitionedGraph=nullptr)
Definition: MeshGraphIO.cpp:53
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:42
std::map< int, TriGeomSharedPtr > TriGeomMap
Definition: TriGeom.h:57
std::map< int, PyrGeomSharedPtr > PyrGeomMap
Definition: PyrGeom.h:76
std::map< int, QuadGeomSharedPtr > QuadGeomMap
Definition: QuadGeom.h:52
std::map< int, SegGeomSharedPtr > SegGeomMap
Definition: SegGeom.h:49
std::map< int, TetGeomSharedPtr > TetGeomMap
Definition: TetGeom.h:84
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:174
std::map< int, PrismGeomSharedPtr > PrismGeomMap
Definition: PrismGeom.h:83
std::shared_ptr< PointGeom > PointGeomSharedPtr
Definition: Geometry.h:57
std::map< int, HexGeomSharedPtr > HexGeomMap
Definition: HexGeom.h:85
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:430
std::map< StdRegions::VarCoeffType, VarCoeffEntry > VarCoeffMap
Definition: StdRegions.hpp:375
StdRegions::ConstFactorMap factors
double NekDouble
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
STL namespace.