File: MPIHelper.cpp

package info (click to toggle)
iqtree 1.5.3%2Bdfsg-2
  • links: PTS, VCS
  • area: main
  • in suites: stretch
  • size: 9,780 kB
  • ctags: 11,529
  • sloc: cpp: 96,162; ansic: 59,874; python: 242; sh: 189; makefile: 45
file content (560 lines) | stat: -rw-r--r-- 17,242 bytes parent folder | download | duplicates (3)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
//
// Created by tung on 6/18/15.
//

#include "MPIHelper.h"
#include "timeutil.h"

/**
 *  Initialize the single getInstance of MPIHelper
 */

MPIHelper& MPIHelper::getInstance() {
    static MPIHelper instance;
#ifndef _IQTREE_MPI
    instance.setProcessID(0);
    instance.setNumProcesses(1);
#endif
    return instance;
}

void MPIHelper::distributeTrees(vector<string> &treeStrings, vector<double> &scores, int tag) {
    if (getNumProcesses() == 1)
        return;
#ifdef _IQTREE_MPI
    vector<int> sourceProcID;
    sourceProcID.insert(sourceProcID.end(), scores.size(), getProcessID());
    TreeCollection outTrees(treeStrings, scores, sourceProcID);
    cleanUpMessages();
    for (int i = 0; i < getNumProcesses(); i++) {
        if (i != getProcessID()) {
            MPI_Request *request = new MPI_Request;
            ObjectStream *os = new ObjectStream(outTrees);
            MPI_Isend(os->getObjectData(), os->getDataLength(), MPI_CHAR, i, tag, MPI_COMM_WORLD, request);
            sentMessages.push_back(make_pair(request, os));
            int flag = 0;
            MPI_Status status;
            MPI_Test(request, &flag, &status);
        }
    }
    //numTreeSent += treeStrings.size();
#endif
}

void MPIHelper::distributeTree(string treeString, double score, int tag) {
    if (getNumProcesses() == 1)
        return;
#ifdef _IQTREE_MPI
    double start = getRealTime();
    vector<string> trees;
    vector<double> scores;
    trees.push_back(treeString);
    scores.push_back(score);
    distributeTrees(trees, scores, tag);
    if (verbose_mode >= VB_MED)
        cout << "Sent tree to other processes in " << getRealTime() - start << " seconds" << endl;
    numTreeSent++;
#endif
}

void MPIHelper::sendTrees(int dest, vector<string> &treeStrings, vector<double> &scores, int tag) {
    if (getNumProcesses() == 1 || dest == getProcessID())
        return;
#ifdef _IQTREE_MPI
    vector<int> sourceProcID;
    sourceProcID.insert(sourceProcID.end(), scores.size(), getProcessID());
    TreeCollection outTrees(treeStrings, scores, sourceProcID);
    cleanUpMessages();
    MPI_Request *request = new MPI_Request;
    ObjectStream *os = new ObjectStream(outTrees);
    MPI_Isend(os->getObjectData(), os->getDataLength(), MPI_CHAR, dest, tag, MPI_COMM_WORLD, request);
    sentMessages.push_back(make_pair(request, os));
    numTreeSent += treeStrings.size();

    int flag = 0;
    MPI_Status status;
    MPI_Test(request, &flag, &status);
#endif
}

void MPIHelper::sendTree(int dest, string treeString, double score, int tag) {
    if (getNumProcesses() == 1 || dest == getProcessID())
        return;
#ifdef _IQTREE_MPI
    StrVector treeStrings;
    treeStrings.push_back(treeString);
    DoubleVector scores;
    scores.push_back(score);
    sendTrees(dest, treeStrings, scores, tag);
#endif
}

int MPIHelper::sendRecvTrees(int dest, vector<string> &treeStrings, vector<double> &scores, int tag) {
    if (getNumProcesses() == 1 || dest == getProcessID())
        return tag;
#ifdef _IQTREE_MPI
    double beginTime = getRealTime();
    // prepare message
    vector<int> sourceProcID;
    sourceProcID.insert(sourceProcID.end(), scores.size(), getProcessID());
    TreeCollection outTrees(treeStrings, scores, sourceProcID);
    ObjectStream *os = new ObjectStream(outTrees);

    // blocking send
    MPI_Send(os->getObjectData(), os->getDataLength(), MPI_CHAR, dest, tag, MPI_COMM_WORLD);
    numTreeSent += treeStrings.size();
    delete os;

    // blocking probe
    MPI_Status status;
    MPI_Probe(dest, MPI_ANY_TAG, MPI_COMM_WORLD, &status);
    int msgCount;
    MPI_Get_count(&status, MPI_CHAR, &msgCount);

    // receive the message
    char *recvBuffer = new char[msgCount];
    MPI_Recv(recvBuffer, msgCount, MPI_CHAR, status.MPI_SOURCE, status.MPI_TAG, MPI_COMM_WORLD, &status);
    treeStrings.clear();
    scores.clear();

    if (status.MPI_TAG != STOP_TAG) {
        os = new ObjectStream(recvBuffer, msgCount);
        TreeCollection curTrees = os->getTreeCollection();
        treeStrings = curTrees.getTreeStrings();
        scores = curTrees.getScores();
        numTreeReceived += treeStrings.size();
    }
    delete [] recvBuffer;

    double endTime = getRealTime();
    cout << "INFO: " << endTime - beginTime << " seconds for " << __func__ << endl;

    return status.MPI_TAG;
#else
    return tag;
#endif
}

int MPIHelper::recvSendTrees(vector<string> &treeStrings, vector<double> &scores, vector<bool> &should_send, int tag) {
    if (getNumProcesses() == 1)
        return 0;
#ifdef _IQTREE_MPI
    double beginTime = getRealTime();
    // blocking probe
    MPI_Status status;
    MPI_Probe(MPI_ANY_SOURCE, MPI_ANY_TAG, MPI_COMM_WORLD, &status);
    int msgCount;
    MPI_Get_count(&status, MPI_CHAR, &msgCount);
    int dest = status.MPI_SOURCE;

    // receive the message
    char *recvBuffer = new char[msgCount];
    MPI_Recv(recvBuffer, msgCount, MPI_CHAR, status.MPI_SOURCE, status.MPI_TAG, MPI_COMM_WORLD, &status);

    // now send message
    if (!should_send[dest]) {
        treeStrings.resize(1, "notree");
        scores.resize(1, -DBL_MAX);
    }
    IntVector sourceProcID;
    sourceProcID.insert(sourceProcID.end(), scores.size(), getProcessID());
    TreeCollection outTrees(treeStrings, scores, sourceProcID);
    ObjectStream *os = new ObjectStream(outTrees);

    // blocking send
    MPI_Send(os->getObjectData(), os->getDataLength(), MPI_CHAR, dest, tag, MPI_COMM_WORLD);
    numTreeSent += treeStrings.size();
    delete os;

    // now extract trees from received buffer
    treeStrings.clear();
    scores.clear();
    os = new ObjectStream(recvBuffer, msgCount);
    TreeCollection curTrees = os->getTreeCollection();
    treeStrings = curTrees.getTreeStrings();
    scores = curTrees.getScores();
    delete [] recvBuffer;
    numTreeReceived += treeStrings.size();
    
    should_send[dest] = false;

    double endTime = getRealTime();
    if (endTime - beginTime > 1)
        cout << "WARNING: " << endTime - beginTime << " seconds for " << __func__ << endl;

    return dest;
#else
    return 0;
#endif
}

void MPIHelper::gatherTrees(TreeCollection &trees) {
    if (getNumProcesses() == 1)
        return;
#ifdef _IQTREE_MPI
    double beginTime = getRealTime();

    if (isMaster()) {
        trees.clear();
        // Master: receive from all Workers
        for (int w = 1; w < getNumProcesses(); w++) {
            // blocking probe
            MPI_Status status;
            MPI_Probe(MPI_ANY_SOURCE, MPI_ANY_TAG, MPI_COMM_WORLD, &status);
            int msgCount;
            MPI_Get_count(&status, MPI_CHAR, &msgCount);
            // receive the message
            char *recvBuffer = new char[msgCount];
            MPI_Recv(recvBuffer, msgCount, MPI_CHAR, status.MPI_SOURCE, status.MPI_TAG, MPI_COMM_WORLD, &status);
            ObjectStream *os = new ObjectStream(recvBuffer, msgCount);
            TreeCollection curTrees = os->getTreeCollection();
            trees.addTrees(curTrees);
            numTreeReceived += curTrees.getNumTrees();
            delete [] recvBuffer;
        }
        cout << trees.getNumTrees() << " trees gathered from workers in ";
    } else {
        // Worker: send trees to Master
        ObjectStream *os = new ObjectStream(trees);
        // blocking send
        MPI_Send(os->getObjectData(), os->getDataLength(), MPI_CHAR, PROC_MASTER, TREE_TAG, MPI_COMM_WORLD);
        numTreeSent += trees.getNumTrees();
        delete os;
        cout << trees.getNumTrees() << " trees sent to master in ";
    }

    double endTime = getRealTime();
    cout << endTime - beginTime << " seconds" << endl;
#endif
}

void MPIHelper::broadcastTrees(TreeCollection &trees) {
    if (getNumProcesses() == 1)
        return;
#ifdef _IQTREE_MPI
    double beginTime = getRealTime();

    // prepare data from Master
    ObjectStream *os;
    int msgCount = 0;
    if (isMaster()) {
        os = new ObjectStream(trees);
        msgCount = os->getDataLength();
    }

    // broadcast the count for workers
    MPI_Bcast(&msgCount, 1, MPI_INT, PROC_MASTER, MPI_COMM_WORLD);

    char *recvBuffer = new char[msgCount];
    if (isMaster())
        memcpy(recvBuffer, os->getObjectData(), msgCount);

    // broadcast trees to workers
    MPI_Bcast(recvBuffer, msgCount, MPI_CHAR, PROC_MASTER, MPI_COMM_WORLD);

    if (isWorker()) {
        os = new ObjectStream(recvBuffer, msgCount);
        trees = os->getTreeCollection();
    }
    delete os;
    delete [] recvBuffer;

    double endTime = getRealTime();
    cout << trees.getNumTrees() << " trees broadcasted to workers in " << endTime - beginTime << " seconds" << endl;

#endif
}


bool MPIHelper::gotMessage() {
    // Check for incoming messages
    if (getNumProcesses() == 1)
        return false;
#ifdef _IQTREE_MPI
    int flag = 0;
    MPI_Status status;
    MPI_Iprobe(MPI_ANY_SOURCE, MPI_ANY_TAG, MPI_COMM_WORLD, &flag, &status);
    if (flag)
        return true;
    else
        return false;
#else
    return false;
#endif
}

void MPIHelper::sendMsg(int tag, string msg) {
    if (getNumProcesses() == 1)
        return;
#ifdef _IQTREE_MPI
    if (tag == STOP_TAG)
        cleanUpMessages();
    for (int i = 0; i < getNumProcesses(); i++) {
        if (i != getProcessID()) {
            MPI_Request *request = new MPI_Request;
            ObjectStream *os = new ObjectStream(msg.c_str(), msg.size()+1);
            MPI_Isend(os->getObjectData(), os->getDataLength(), MPI_CHAR, i, tag, MPI_COMM_WORLD, request);
            sentMessages.push_back(make_pair(request, os));
            int flag = 0;
            MPI_Status status;
            MPI_Test(request, &flag, &status);
        }
    }
#endif
}

bool MPIHelper::checkMsg(int tag, string &msg) {
    if (getNumProcesses() == 1)
        return true;
#ifdef _IQTREE_MPI
    int flag=0;
    MPI_Status status;
    char *recvBuffer;
    int numBytes;
    // Check for incoming messages
    MPI_Iprobe(PROC_MASTER, tag, MPI_COMM_WORLD, &flag, &status);
    // flag == true if there is a message
    if (flag) {
        MPI_Get_count(&status, MPI_CHAR, &numBytes);
        recvBuffer = new char[numBytes];
        MPI_Recv(recvBuffer, numBytes, MPI_CHAR, status.MPI_SOURCE, status.MPI_TAG, MPI_COMM_WORLD, &status);
        msg = recvBuffer;
        delete[] recvBuffer;
        return true;
    }
#endif
    return false;
}

bool MPIHelper::checkMsg(int tag) {
    if (getNumProcesses() == 1) {
        return false;
    }
#ifdef _IQTREE_MPI
    string msg;
    if (checkMsg(tag, msg)) {
        cout << "Worker " << getProcessID() << " gets message " << msg << endl;
        return true;
    }
#endif
    return false;
}


void MPIHelper::receiveTrees(bool fromAll, int maxNumTrees, TreeCollection &trees, int tag) {
    if (getNumProcesses() == 1) {
        return;
    }
#ifdef _IQTREE_MPI
    int flag = 0;
    int minNumTrees = 0;
    bool nodes[getNumProcesses()];
    if (fromAll)
        minNumTrees = getNumProcesses() - 1;
    for (int i = 0; i < getNumProcesses(); i++)
        nodes[i] = false;
    nodes[getProcessID()] = true;
    // Process all pending messages
    MPI_Status status;
    size_t totalMsgSize = 0;
    do {
        char* recvBuffer;
        int numBytes;
        flag = 0;
        // Check for incoming messages
        MPI_Iprobe(MPI_ANY_SOURCE, tag, MPI_COMM_WORLD, &flag, &status);
        // flag == true if there is a message
        if (flag) {
            //cout << "Getting messages from node " << status.MPI_SOURCE << endl;
            MPI_Get_count(&status, MPI_CHAR, &numBytes);
            totalMsgSize += numBytes;
            recvBuffer = new char[numBytes];
            MPI_Recv(recvBuffer, numBytes, MPI_CHAR, status.MPI_SOURCE, status.MPI_TAG, MPI_COMM_WORLD, &status);
            ObjectStream os(recvBuffer, numBytes);
            if (status.MPI_TAG == STOP_TAG) {
                cout <<  os.getObjectData() << endl;
                MPI_Finalize();
                exit(0);
            }
            TreeCollection curTrees = os.getTreeCollection();
            trees.addTrees(curTrees);
            if (trees.getNumTrees() >= maxNumTrees) {
                break;
            }
            if (fromAll && !nodes[status.MPI_SOURCE]) {
                nodes[status.MPI_SOURCE] = true;
                minNumTrees--;
            }
            delete [] recvBuffer;
        }
    } while (minNumTrees > 0 || flag);
    numTreeReceived += trees.getNumTrees();
    if (trees.getNumTrees() > 0) {
        cout << "Proc " << getProcessID() << ": " << trees.getNumTrees() << " trees received from other processes (" << totalMsgSize << " bytes)" << endl;
    }
#endif
}

int MPIHelper::receiveTrees(TreeCollection &trees, int tag) {
    if (getNumProcesses() == 1) {
        return -1;
    }
#ifdef _IQTREE_MPI
    int flag = 0;
    // Process all pending messages
    MPI_Status status;
    char* recvBuffer;
    int numBytes;
    // Check for incoming messages
    MPI_Iprobe(MPI_ANY_SOURCE, tag, MPI_COMM_WORLD, &flag, &status);
    // flag == true if there is a message
    if (flag) {
        //cout << "Getting messages from node " << status.MPI_SOURCE << endl;
        MPI_Get_count(&status, MPI_CHAR, &numBytes);
        recvBuffer = new char[numBytes];
        MPI_Recv(recvBuffer, numBytes, MPI_CHAR, status.MPI_SOURCE, status.MPI_TAG, MPI_COMM_WORLD, &status);
        ObjectStream os(recvBuffer, numBytes);
        TreeCollection curTrees = os.getTreeCollection();
        trees.addTrees(curTrees);
        delete [] recvBuffer;
        return status.MPI_SOURCE;
    }
#endif
    return -1;
}

int MPIHelper::cleanUpMessages() {
#ifdef _IQTREE_MPI
    int numMsgCleaned = 0;
    // change iterator to index to avoid iterator being invalidated after erase()
    for (int i = 0; i < sentMessages.size(); ) {
        int flag = 0;
        MPI_Status status;
        MPI_Test(sentMessages[i].first, &flag, &status);
        if (flag) {
            delete sentMessages[i].first;
            delete sentMessages[i].second;
            numMsgCleaned++;
            sentMessages.erase(sentMessages.begin()+i);
        } else {
            i++;
        }
    }
    if (verbose_mode >= VB_MED && numMsgCleaned)
        cout << numMsgCleaned << " messages sent and cleaned up" << endl;
    return numMsgCleaned;
#else
    return 0;
#endif
}

#ifdef _IQTREE_MPI
void MPIHelper::sendString(string &str, int dest, int tag) {
    char *buf = (char*)str.c_str();
    MPI_Send(buf, str.length()+1, MPI_CHAR, dest, tag, MPI_COMM_WORLD);
}

void MPIHelper::sendCheckpoint(Checkpoint *ckp, int dest) {
    stringstream ss;
    ckp->dump(ss);
    string str = ss.str();
    sendString(str, dest, TREE_TAG);
}


int MPIHelper::recvString(string &str, int src, int tag) {
    MPI_Status status;
    MPI_Probe(src, tag, MPI_COMM_WORLD, &status);
    int msgCount;
    MPI_Get_count(&status, MPI_CHAR, &msgCount);
    // receive the message
    char *recvBuffer = new char[msgCount];
    MPI_Recv(recvBuffer, msgCount, MPI_CHAR, status.MPI_SOURCE, status.MPI_TAG, MPI_COMM_WORLD, &status);
    str = recvBuffer;
    delete [] recvBuffer;
    return status.MPI_SOURCE;
}

int MPIHelper::recvCheckpoint(Checkpoint *ckp, int src) {
    string str;
    int proc = recvString(str, src, TREE_TAG);
    stringstream ss(str);
    ckp->load(ss);
    return proc;
}

void MPIHelper::broadcastCheckpoint(Checkpoint *ckp) {
    int msgCount = 0;
    stringstream ss;
    string str;
    if (isMaster()) {
        ckp->dump(ss);
        str = ss.str();
        msgCount = str.length()+1;
    }

    // broadcast the count for workers
    MPI_Bcast(&msgCount, 1, MPI_INT, PROC_MASTER, MPI_COMM_WORLD);

    char *recvBuffer = new char[msgCount];
    if (isMaster())
        memcpy(recvBuffer, str.c_str(), msgCount);

    // broadcast trees to workers
    MPI_Bcast(recvBuffer, msgCount, MPI_CHAR, PROC_MASTER, MPI_COMM_WORLD);

    if (isWorker()) {
        ss.clear();
        ss.str(recvBuffer);
        ckp->load(ss);
    }
    delete [] recvBuffer;
}

void MPIHelper::gatherCheckpoint(Checkpoint *ckp) {
    stringstream ss;
    ckp->dump(ss);
    string str = ss.str();
    int msgCount = str.length();

    // first send the counts to MASTER
    int *msgCounts = NULL, *displ = NULL;
    char *recvBuffer = NULL;
    int totalCount = 0;

    if (isMaster()) {
        msgCounts = new int[getNumProcesses()];
        displ = new int[getNumProcesses()];
    }
    MPI_Gather(&msgCount, 1, MPI_INT, msgCounts, 1, MPI_INT, PROC_MASTER, MPI_COMM_WORLD);

    // now real contents to MASTER
    if (isMaster()) {
        for (int i = 0; i < getNumProcesses(); i++) {
            displ[i] = totalCount;
            totalCount += msgCounts[i];
        }
        recvBuffer = new char[totalCount+1];
        memset(recvBuffer, 0, totalCount+1);
    }
    char *buf = (char*)str.c_str();
    MPI_Gatherv(buf, msgCount, MPI_CHAR, recvBuffer, msgCounts, displ, MPI_CHAR, PROC_MASTER, MPI_COMM_WORLD);

    if (isMaster()) {
        // now decode the buffer
        ss.clear();
        ss.str(recvBuffer);
        ckp->load(ss);

        delete [] recvBuffer;
        delete [] displ;
        delete [] msgCounts;
    }
}

#endif

MPIHelper::~MPIHelper() {
//    cleanUpMessages();
}