gmsh-TingyuanDoc  0.1
An Open-Source Timing-driven Analytical Mixed-size FPGA Placer
GmshRemote.cpp
Go to the documentation of this file.
1 // Gmsh - Copyright (C) 1997-2022 C. Geuzaine, J.-F. Remacle
2 //
3 // See the LICENSE.txt file in the Gmsh root directory for license information.
4 // Please report all issues on https://gitlab.onelab.info/gmsh/gmsh/issues.
5 
6 #include "GmshConfig.h"
7 
8 #if defined(HAVE_MPI)
9 #include <mpi.h>
10 #define MPI_GMSH_COMPUTE_VIEW 1
11 #define MPI_GMSH_DATA_READY 2
12 #define MPI_GMSH_VARRAY 3
13 #define MPI_GMSH_VARRAY_LEN 4
14 #define MPI_GMSH_SHUTDOWN 5
15 #define MPI_GMSH_PARSE_STRING 6
16 #define MPI_GMSH_MERGE_FILE 7
17 #endif
18 
19 #include <sstream>
20 #include "GmshMessage.h"
21 
22 #if defined(HAVE_ONELAB) && defined(HAVE_POST)
23 
24 #include "onelab.h"
25 #include "OpenFile.h"
26 #include "OS.h"
27 #include "VertexArray.h"
28 #include "GmshRemote.h"
29 #include "PView.h"
30 #include "PViewOptions.h"
31 #include "PViewData.h"
32 #include "PViewDataRemote.h"
33 
34 static void ComputeAndSendVertexArrays(GmshClient *client, bool compute = true)
35 {
36  for(std::size_t i = 0; i < PView::list.size(); i++) {
37  PView *p = PView::list[i];
38  if(compute) p->fillVertexArrays();
39  PViewData *data = p->getData();
40  PViewOptions *opt = p->getOptions();
41  double min = data->getMin(), max = data->getMax();
43  min = data->getMin(opt->timeStep);
44  max = data->getMax(opt->timeStep);
45  }
46  VertexArray *va[4] = {p->va_points, p->va_lines, p->va_triangles,
47  p->va_vectors};
48  for(int type = 0; type < 4; type++) {
49  if(va[type]) {
50  int len;
51  char *str = va[type]->toChar(p->getTag(), data->getName(), type + 1,
52  min, max, data->getNumTimeSteps(),
53  data->getTime(opt->timeStep),
54  data->getBoundingBox(), len);
55  client->SendMessage(GmshSocket::GMSH_VERTEX_ARRAY, len, str);
56  delete[] str;
57  }
58  }
59  }
60 }
61 
62 #if defined(HAVE_MPI)
63 // This version sends VArrays using MPI
64 static void ComputeAndSendVertexArrays()
65 {
66  // compute...
67  for(std::size_t i = 0; i < PView::list.size(); i++)
68  PView::list[i]->fillVertexArrays();
69 
70  // ...and send
71  int nbArrays = PView::list.size() * 4;
72  MPI_Send(&nbArrays, 1, MPI_INT, 0, MPI_GMSH_DATA_READY, MPI_COMM_WORLD);
73 
74  for(std::size_t i = 0; i < PView::list.size(); i++) {
75  PView *p = PView::list[i];
76  PViewData *data = p->getData();
77  PViewOptions *opt = p->getOptions();
78  double min = data->getMin(), max = data->getMax();
80  min = data->getMin(opt->timeStep);
81  max = data->getMax(opt->timeStep);
82  }
83  VertexArray *va[4] = {p->va_points, p->va_lines, p->va_triangles,
84  p->va_vectors};
85  for(int type = 0; type < 4; type++) {
86  if(va[type]) {
87  int len;
88  char *str = va[type]->toChar(p->getTag(), data->getName(), type + 1,
89  min, max, data->getNumTimeSteps(),
90  data->getTime(opt->timeStep),
91  data->getBoundingBox(), len);
92  MPI_Send(&len, 1, MPI_INT, 0, MPI_GMSH_VARRAY_LEN, MPI_COMM_WORLD);
93  MPI_Send(str, len, MPI_CHAR, 0, MPI_GMSH_VARRAY, MPI_COMM_WORLD);
94  delete[] str;
95  }
96  }
97  }
98 }
99 
100 // Merge the vertex arrays
101 static void AddToVertexArrays(int length, const char *bytes, int swap)
102 {
103  std::string name;
104  int num, type, numSteps;
105  double min, max, time, xmin, ymin, zmin, xmax, ymax, zmax;
106  VertexArray::decodeHeader(length, bytes, swap, name, num, type, min, max,
107  numSteps, time, xmin, ymin, zmin, xmax, ymax, zmax);
108 
109  PView *p = PView::list[num - 1];
110  PViewData *data = p->getData();
111 
112  VertexArray *varrays[4] = {p->va_points, p->va_lines, p->va_triangles,
113  p->va_vectors};
114 
115  VertexArray *va = varrays[type - 1];
116 
117  if(data->getMin() > min) data->setMin(min);
118  if(data->getMax() < max) data->setMax(max);
119 
120  SBoundingBox3d bbox(xmin, ymin, zmin, xmax, ymax, zmax);
121  SBoundingBox3d bb = data->getBoundingBox();
122  bb += bbox;
123 
124  data->setBoundingBox(bb);
125 
126  if(type == 4) type = 2;
127  VertexArray *toAdd = new VertexArray(type, 100);
128  toAdd->fromChar(length, bytes, swap);
129  va->merge(toAdd);
130  delete toAdd;
131 }
132 #endif
133 
134 static void GatherAndSendVertexArrays(GmshClient *client, bool swap)
135 {
136 #if defined(HAVE_MPI)
137  // int rank = Msg::GetCommRank();
138  int nbDaemon = Msg::GetCommSize();
139  // tell every node to start computing
140  int mpi_msg = MPI_GMSH_COMPUTE_VIEW;
141  MPI_Bcast(&mpi_msg, 1, MPI_INT, 0, MPI_COMM_WORLD);
142  // fill the arrays on the master node
143  for(std::size_t i = 0; i < PView::list.size(); i++)
144  PView::list[i]->fillVertexArrays();
145  // wait and send the data from every other node
146  for(int i = 0; i < nbDaemon - 1; i++) {
147  int nbArrays;
148  MPI_Status status;
149  MPI_Recv(&nbArrays, 1, MPI_INT, MPI_ANY_SOURCE, MPI_GMSH_DATA_READY,
150  MPI_COMM_WORLD, &status);
151  // int source = status.MPI_SOURCE;
152  // get each varray in turn, then add it to the varrays of
153  // the master node
154  for(int j = 0; j < nbArrays; j++) {
155  int len;
156  MPI_Status status2;
157  MPI_Recv(&len, 1, MPI_INT, status.MPI_SOURCE, MPI_GMSH_VARRAY_LEN,
158  MPI_COMM_WORLD, &status2);
159  char str[len];
160  MPI_Recv(str, len, MPI_CHAR, status.MPI_SOURCE, MPI_GMSH_VARRAY,
161  MPI_COMM_WORLD, &status2);
162  AddToVertexArrays(len, str, swap);
163  }
164  }
165  ComputeAndSendVertexArrays(client, false);
166 #endif
167 }
168 
169 int GmshRemote()
170 {
171  GmshClient *client = Msg::GetGmshClient();
172 
173  int rank = Msg::GetCommRank();
174  int nbDaemon = Msg::GetCommSize();
175 
176  if(!client && rank == 0) return 0;
177 
178  if(client && nbDaemon < 2)
179  ComputeAndSendVertexArrays(client);
180  else if(client && nbDaemon >= 2 && rank == 0)
181  GatherAndSendVertexArrays(client, false);
182 
183  while(1) {
184  // on the node with MPI rank 0, communicate through a socket
185  if(rank == 0) {
186  // stop if we have no communications for 5 minutes
187  int ret = client->Select(300, 0);
188  if(!ret) {
189  client->Info("Timeout: stopping remote Gmsh...");
190  break;
191  }
192  else if(ret < 0) {
193  client->Error("Error on select: stopping remote Gmsh...");
194  break;
195  }
196 
197  int type, length, swap;
198  if(!client->ReceiveHeader(&type, &length, &swap)) {
199  client->Error(
200  "Did not receive message header: stopping remote Gmsh...");
201  break;
202  }
203 
204  char *msg = new char[length + 1];
205  if(!client->ReceiveString(length, msg)) {
206  client->Error("Did not receive message body: stopping remote Gmsh...");
207  delete[] msg;
208  break;
209  }
210 
211  if(type == GmshSocket::GMSH_STOP) {
212  client->Info("Stopping remote Gmsh...");
213  delete[] msg;
214  break;
215  }
216  else if(type == GmshSocket::GMSH_VERTEX_ARRAY) {
217  ParseString(msg);
218 #if !defined(HAVE_MPI)
219  ComputeAndSendVertexArrays(client);
220 #else
221  int mpi_msg = MPI_GMSH_PARSE_STRING;
222  MPI_Bcast(&mpi_msg, 1, MPI_INT, 0, MPI_COMM_WORLD);
223  MPI_Bcast(&length, 1, MPI_INT, 0, MPI_COMM_WORLD);
224  MPI_Bcast(msg, length, MPI_CHAR, 0, MPI_COMM_WORLD);
225  GatherAndSendVertexArrays(client, swap);
226 #endif
227  }
228  else if(type == GmshSocket::GMSH_MERGE_FILE) {
229  MergeFile(msg);
230 #if !defined(HAVE_MPI)
231  ComputeAndSendVertexArrays(client);
232 #else
233  int mpi_msg = MPI_GMSH_MERGE_FILE;
234  MPI_Bcast(&mpi_msg, 1, MPI_INT, 0, MPI_COMM_WORLD);
235  MPI_Bcast(&length, 1, MPI_INT, 0, MPI_COMM_WORLD);
236  MPI_Bcast(msg, length, MPI_CHAR, 0, MPI_COMM_WORLD);
237  GatherAndSendVertexArrays(client, swap);
238 #endif
239  }
240  else if(type == GmshSocket::GMSH_PARSE_STRING) {
241  ParseString(msg);
242 #if defined(HAVE_MPI)
243  int mpi_msg = MPI_GMSH_PARSE_STRING;
244  MPI_Bcast(&mpi_msg, 1, MPI_INT, 0, MPI_COMM_WORLD);
245  MPI_Bcast(&length, 1, MPI_INT, 0, MPI_COMM_WORLD);
246  MPI_Bcast(msg, length, MPI_CHAR, 0, MPI_COMM_WORLD);
247 #endif
248  }
249  else if(type == GmshSocket::GMSH_SPEED_TEST) {
250  client->Info("Sending huge array");
251  std::string huge(500000000, 'a');
252  client->SpeedTest(huge.c_str());
253  }
254  else {
255  client->Error("Ignoring unknown message");
256  }
257 
258  delete[] msg;
259  }
260  else { // if we're not on the master node (rank != 0) wait for him...
261 #if defined(HAVE_MPI)
262  int mpi_msg;
263  MPI_Bcast(&mpi_msg, 1, MPI_INT, 0, MPI_COMM_WORLD);
264  if(mpi_msg == MPI_GMSH_COMPUTE_VIEW)
265  ComputeAndSendVertexArrays();
266  else if(mpi_msg == MPI_GMSH_SHUTDOWN)
267  Msg::Exit(0);
268  else if(mpi_msg == MPI_GMSH_PARSE_STRING) {
269  int length;
270  MPI_Bcast(&length, 1, MPI_INT, 0, MPI_COMM_WORLD);
271  char msg[length];
272  MPI_Bcast(msg, length, MPI_CHAR, 0, MPI_COMM_WORLD);
273  ParseString(msg);
274  }
275  else if(mpi_msg == MPI_GMSH_MERGE_FILE) {
276  int length;
277  MPI_Bcast(&length, 1, MPI_INT, 0, MPI_COMM_WORLD);
278  char msg[length];
279  MPI_Bcast(msg, length, MPI_CHAR, 0, MPI_COMM_WORLD);
280  MergeFile(msg);
281  }
282 #endif
283  }
284  }
285 
286 #if defined(HAVE_MPI)
287  int mpi_msg = MPI_GMSH_SHUTDOWN;
288  MPI_Bcast(&mpi_msg, 1, MPI_INT, 0, MPI_COMM_WORLD);
289 #endif
290 
291  return 0;
292 }
293 
294 #else
295 
297 {
298  Msg::Error("GmshRemote requires Post and ONELAB modules");
299  return 0;
300 }
301 
302 #endif
GmshSocket::GMSH_VERTEX_ARRAY
@ GMSH_VERTEX_ARRAY
Definition: GmshSocket.h:74
PViewOptions::timeStep
int timeStep
Definition: PViewOptions.h:61
PView::va_points
VertexArray * va_points
Definition: PView.h:156
GmshSocket::SendMessage
void SendMessage(int type, int length, const void *msg)
Definition: GmshSocket.h:186
PView
Definition: PView.h:27
VertexArray::merge
void merge(VertexArray *va)
Definition: VertexArray.cpp:343
GmshSocket::ReceiveString
int ReceiveString(int len, char *str)
Definition: GmshSocket.h:235
VertexArray::fromChar
void fromChar(int length, const char *bytes, int swap)
Definition: VertexArray.cpp:313
GmshSocket::ReceiveHeader
int ReceiveHeader(int *type, int *len, int *swap)
Definition: GmshSocket.h:212
MergeFile
int MergeFile(const std::string &fileName, bool errorIfMissing, bool setBoundingBox, bool importPhysicalsInOnelab, int partitionToRead)
Definition: OpenFile.cpp:298
GmshRemote
int GmshRemote()
Definition: GmshRemote.cpp:296
GmshSocket::SpeedTest
void SpeedTest(const char *str)
Definition: GmshSocket.h:205
VertexArray::toChar
char * toChar(int num, const std::string &name, int type, double min, double max, int numsteps, double time, const SBoundingBox3d &bbox, int &len)
Definition: VertexArray.cpp:234
PViewData::getNumTimeSteps
virtual int getNumTimeSteps()=0
OS.h
GmshSocket::GMSH_PARSE_STRING
@ GMSH_PARSE_STRING
Definition: GmshSocket.h:73
Msg::Error
static void Error(const char *fmt,...)
Definition: GmshMessage.cpp:482
VertexArray.h
VertexArray
Definition: VertexArray.h:151
PView.h
Msg::GetGmshClient
static GmshClient * GetGmshClient()
Definition: GmshMessage.cpp:1272
PView::va_vectors
VertexArray * va_vectors
Definition: PView.h:156
GmshSocket::GMSH_SPEED_TEST
@ GMSH_SPEED_TEST
Definition: GmshSocket.h:82
GmshMessage.h
PViewData.h
VertexArray::decodeHeader
static int decodeHeader(int length, const char *bytes, int swap, std::string &name, int &tag, int &type, double &min, double &max, int &numSteps, double &time, double &xmin, double &ymin, double &zmin, double &xmax, double &ymax, double &zmax)
Definition: VertexArray.cpp:273
GmshClient
Definition: GmshSocket.h:261
ParseString
void ParseString(const std::string &str, bool inCurrentModelDir)
Definition: OpenFile.cpp:259
PViewData::getTime
virtual double getTime(int step)
Definition: PViewData.h:94
PViewData::getMax
virtual double getMax(int step=-1, bool onlyVisible=false, int tensorRep=0, int forceNumComponents=0, int componentMap[9]=nullptr)=0
PView::getTag
int getTag()
Definition: PView.h:89
Msg::GetCommRank
static int GetCommRank()
Definition: GmshMessage.cpp:219
GmshRemote.h
PViewData::getBoundingBox
virtual SBoundingBox3d getBoundingBox(int step=-1)=0
swap
void swap(double &a, double &b)
Definition: meshTriangulation.cpp:27
PViewData::setMin
virtual void setMin(double min)=0
PView::getData
PViewData * getData(bool useAdaptiveIfAvailable=false)
Definition: PView.cpp:233
PViewOptions::rangeType
int rangeType
Definition: PViewOptions.h:59
onelab.h
PViewOptions.h
GmshSocket::Info
void Info(const char *str)
Definition: GmshSocket.h:198
Msg::Exit
static void Exit(int level)
Definition: GmshMessage.cpp:384
PViewDataRemote.h
GmshSocket::GMSH_STOP
@ GMSH_STOP
Definition: GmshSocket.h:67
Msg::GetCommSize
static int GetCommSize()
Definition: GmshMessage.cpp:224
PView::fillVertexArrays
bool fillVertexArrays()
Definition: PViewVertexArrays.cpp:1568
GmshSocket::Select
int Select(int seconds, int microseconds, int socket=-1)
Definition: GmshSocket.h:173
PViewData
Definition: PViewData.h:29
length
double length(Quaternion &q)
Definition: Camera.cpp:346
PViewOptions::PerTimeStep
@ PerTimeStep
Definition: PViewOptions.h:37
PView::va_triangles
VertexArray * va_triangles
Definition: PView.h:156
PViewData::setMax
virtual void setMax(double max)=0
PViewOptions
Definition: PViewOptions.h:16
PView::getOptions
PViewOptions * getOptions()
Definition: PView.h:81
PViewData::getMin
virtual double getMin(int step=-1, bool onlyVisible=false, int tensorRep=0, int forceNumComponents=0, int componentMap[9]=nullptr)=0
PViewData::getName
virtual std::string getName()
Definition: PViewData.h:70
GmshSocket::GMSH_MERGE_FILE
@ GMSH_MERGE_FILE
Definition: GmshSocket.h:72
PViewData::setBoundingBox
virtual void setBoundingBox(SBoundingBox3d &box)=0
PView::list
static std::vector< PView * > list
Definition: PView.h:112
PView::va_lines
VertexArray * va_lines
Definition: PView.h:156
SBoundingBox3d
Definition: SBoundingBox3d.h:21
GmshSocket::Error
void Error(const char *str)
Definition: GmshSocket.h:200
OpenFile.h