gmsh-TingyuanDoc  0.1
An Open-Source Timing-driven Analytical Mixed-size FPGA Placer
dofManager.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 #ifdef HAVE_MPI
9 #include "mpi.h"
10 #endif
11 
12 #include <dofManager.h>
13 
14 template <> void dofManager<double>::scatterSolution()
15 {
16 #ifdef HAVE_MPI
17  if(!_parallelFinalized) {
19  }
20  MPI_Status status;
21  std::vector<std::vector<double> > sendBuf(Msg::GetCommSize()),
22  recvBuf(Msg::GetCommSize());
23  std::vector<MPI_Request> reqRecv(Msg::GetCommSize()),
24  reqSend(Msg::GetCommSize());
25  for(int i = 0; i < Msg::GetCommSize(); i++) {
26  reqRecv[i] = reqSend[i] = MPI_REQUEST_NULL;
27  if(!ghostByProc[i].empty()) {
28  recvBuf[i].resize(ghostByProc[i].size());
29  MPI_Irecv(&recvBuf[i][0], recvBuf[i].size(), MPI_DOUBLE, i, 0,
30  MPI_COMM_WORLD, &reqRecv[i]);
31  }
32  if(!parentByProc[i].empty()) {
33  getDofValue(parentByProc[i], sendBuf[i]);
34  MPI_Isend(&sendBuf[i][0], sendBuf[i].size(), MPI_DOUBLE, i, 0,
35  MPI_COMM_WORLD, &reqSend[i]);
36  }
37  }
38  int index;
39  while(MPI_Waitany(Msg::GetCommSize(), &reqRecv[0], &index, &status) == 0 &&
40  index != MPI_UNDEFINED) {
41  if(status.MPI_TAG == 0)
42  for(std::size_t j = 0; j < recvBuf[index].size(); j++) {
43  ghostValue[ghostByProc[index][j]] = recvBuf[index][j];
44  // const Dof &dof = ghostByProc[index][j];
45  }
46  }
47  MPI_Waitall(Msg::GetCommSize(), &reqSend[0], MPI_STATUS_IGNORE);
48 #endif
49 }
50 
52 {
53  _localSize = unknown.size();
54 #ifdef HAVE_MPI
55  int numStart;
56  int numTotal;
57  MPI_Status status;
59  ghostByProc.resize(Msg::GetCommSize());
60  if(Msg::GetCommRank() == 0) {
61  numStart = 0;
62  }
63  else
64  MPI_Recv(&numStart, 1, MPI_INT, Msg::GetCommRank() - 1, 0, MPI_COMM_WORLD,
65  &status);
66  numTotal = numStart + _localSize;
67  if(Msg::GetCommRank() != Msg::GetCommSize() - 1)
68  MPI_Send(&numTotal, 1, MPI_INT, Msg::GetCommRank() + 1, 0, MPI_COMM_WORLD);
69  MPI_Bcast(&numTotal, 1, MPI_INT, Msg::GetCommSize() - 1, MPI_COMM_WORLD);
70  for(auto it = unknown.begin(); it != unknown.end(); it++)
71  it->second += numStart;
72  std::vector<std::list<Dof> > ghostedByProc;
73  int *nRequest = new int[Msg::GetCommSize()];
74  int *nRequested = new int[Msg::GetCommSize()];
75  for(int i = 0; i < Msg::GetCommSize(); i++) nRequest[i] = 0;
76  for(auto it = ghostByDof.begin(); it != ghostByDof.end(); it++) {
77  int procId = it->second.first;
78  it->second.second = nRequest[procId]++;
79  }
80  MPI_Alltoall(nRequest, 1, MPI_INT, nRequested, 1, MPI_INT, MPI_COMM_WORLD);
81  long int **recv0 = new long int *[Msg::GetCommSize()];
82  int **recv1 = new int *[Msg::GetCommSize()];
83  long int **send0 = new long int *[Msg::GetCommSize()];
84  int **send1 = new int *[Msg::GetCommSize()];
85  MPI_Request *reqRecv0 = new MPI_Request[2 * Msg::GetCommSize()];
86  MPI_Request *reqRecv1 = reqRecv0 + Msg::GetCommSize();
87  MPI_Request *reqSend0 = new MPI_Request[Msg::GetCommSize()];
88  MPI_Request *reqSend1 = new MPI_Request[Msg::GetCommSize()];
89  for(int i = 0; i < Msg::GetCommSize(); i++) {
90  send0[i] = new long int[nRequest[i] * 2];
91  recv0[i] = new long int[nRequested[i] * 2];
92  send1[i] = new int[nRequested[i]];
93  recv1[i] = new int[nRequest[i]];
94  reqSend0[i] = reqSend1[i] = reqRecv0[i] = reqRecv1[i] = MPI_REQUEST_NULL;
95  parentByProc[i].resize(nRequested[i], Dof(0, 0));
96  ghostByProc[i].resize(nRequest[i], Dof(0, 0));
97  }
98  for(int i = 0; i < Msg::GetCommSize(); i++) nRequest[i] = 0;
99  for(auto it = ghostByDof.begin(); it != ghostByDof.end(); it++) {
100  int proc = it->second.first;
101  send0[proc][nRequest[proc] * 2] = it->first.getEntity();
102  send0[proc][nRequest[proc] * 2 + 1] = it->first.getType();
103  ghostByProc[proc][nRequest[proc]] = it->first;
104  nRequest[proc]++;
105  }
106  for(int i = 0; i < Msg::GetCommSize(); i++) {
107  if(nRequested[i] > 0) {
108  MPI_Irecv(recv0[i], 2 * nRequested[i], MPI_LONG, i, 0, MPI_COMM_WORLD,
109  &reqRecv0[i]);
110  }
111  if(nRequest[i] > 0) {
112  MPI_Irecv(recv1[i], 2 * nRequest[i], MPI_INT, i, 1, MPI_COMM_WORLD,
113  &reqRecv1[i]);
114  MPI_Isend(send0[i], 2 * nRequest[i], MPI_LONG, i, 0, MPI_COMM_WORLD,
115  &reqSend0[i]);
116  }
117  }
118  int index;
119  while(MPI_Waitany(2 * Msg::GetCommSize(), reqRecv0, &index, &status) == 0 &&
120  index != MPI_UNDEFINED) {
121  if(status.MPI_TAG == 0) {
122  for(int j = 0; j < nRequested[index]; j++) {
123  Dof d(recv0[index][j * 2], recv0[index][j * 2 + 1]);
124  auto it = unknown.find(d);
125  if(it == unknown.end())
126  Msg::Error("ghost Dof does not exist on parent process");
127  send1[index][j] = it->second;
128  parentByProc[index][j] = d;
129  }
130  MPI_Isend(send1[index], nRequested[index], MPI_INT, index, 1,
131  MPI_COMM_WORLD, &reqSend1[index]);
132  }
133  }
134  for(int i = 0; i < Msg::GetCommSize(); i++)
135  for(int i = 0; i < Msg::GetCommSize(); i++) nRequest[i] = 0;
136  for(auto it = ghostByDof.begin(); it != ghostByDof.end(); it++) {
137  int proc = it->second.first;
138  unknown[it->first] = recv1[proc][nRequest[proc]++];
139  }
140  MPI_Waitall(Msg::GetCommSize(), reqSend0, MPI_STATUS_IGNORE);
141  MPI_Waitall(Msg::GetCommSize(), reqSend1, MPI_STATUS_IGNORE);
142  for(int i = 0; i < Msg::GetCommSize(); i++) {
143  delete[] send0[i];
144  delete[] send1[i];
145  delete[] recv0[i];
146  delete[] recv1[i];
147  }
148  delete[] send0;
149  delete[] send1;
150  delete[] recv0;
151  delete[] recv1;
152  delete[] reqSend0;
153  delete[] reqSend1;
154  delete[] reqRecv0;
155 #endif
156  _parallelFinalized = true;
157 }
dofManagerBase::unknown
std::map< Dof, int > unknown
Definition: dofManager.h:89
Msg::Error
static void Error(const char *fmt,...)
Definition: GmshMessage.cpp:482
dofManager::ghostValue
std::map< Dof, T > ghostValue
Definition: dofManager.h:135
dofManagerBase::parentByProc
std::vector< std::vector< Dof > > parentByProc
Definition: dofManager.h:98
dofManager::scatterSolution
void scatterSolution()
dofManagerBase::_parallelFinalize
void _parallelFinalize()
Definition: dofManager.cpp:51
dofManagerBase::_localSize
int _localSize
Definition: dofManager.h:99
Dof
Definition: dofManager.h:19
Msg::GetCommRank
static int GetCommRank()
Definition: GmshMessage.cpp:219
dofManager.h
Msg::GetCommSize
static int GetCommSize()
Definition: GmshMessage.cpp:224
dofManagerBase::_parallelFinalized
bool _parallelFinalized
Definition: dofManager.h:100
dofManager::getDofValue
virtual void getDofValue(std::vector< Dof > &keys, std::vector< dataVec > &Vals)
Definition: dofManager.h:235
dofManagerBase::ghostByProc
std::vector< std::vector< Dof > > ghostByProc
Definition: dofManager.h:98
dofManagerBase::ghostByDof
std::map< Dof, std::pair< int, int > > ghostByDof
Definition: dofManager.h:97