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