File: Wave2d.cpp

package info (click to toggle)
mpich 4.3.2-2
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid
  • size: 101,184 kB
  • sloc: ansic: 1,040,629; cpp: 82,270; javascript: 40,763; perl: 27,933; python: 16,041; sh: 14,676; xml: 14,418; f90: 12,916; makefile: 9,270; fortran: 8,046; java: 4,635; asm: 324; ruby: 103; awk: 27; lisp: 19; php: 8; sed: 4
file content (191 lines) | stat: -rw-r--r-- 6,675 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
#include "mpi.h"
#include <cmath>
#include <cassert>
#include <cstdlib>
#include <cstdio>
#include <ctime>

#define LIN(x, y, dimX) ((x) + (y)*(dimX))
#define LIN2(x,y,dimX) (LIN(x+1,y+1,dimX+2))
#define CYCLIN(x, y, dimX, dimY) LIN((x+dimX) % (dimX), (y+dimY) % (dimY), (dimX))

static double ks = 1.0;

class Grid {
  private:
    double *myGridPrev_, *myGridCurr_, *myGridNext_;
    double *leftEdgeIn_, *leftEdgeOut_, *rightEdgeIn_, *rightEdgeOut_;
    int numGridsX_, numGridsY_, dimX_, dimY_;
    int myGridX_, myGridY_;

  public:
    Grid(int dimX, int dimY, int numGridsX, int numGridsY, int numInitialPerturbations)
        : dimX_(dimX), dimY_(dimY), numGridsX_(numGridsX), numGridsY_(numGridsY) {
        myGridPrev_ = new double[(dimX+2)*(2+dimY)]();
        myGridCurr_ = new double[(dimX+2)*(2+dimY)]();
        myGridNext_ = new double[(dimX+2)*(2+dimY)]();
        leftEdgeIn_   = new double[dimY];
        rightEdgeIn_  = new double[dimY];
        leftEdgeOut_  = new double[dimY];
        rightEdgeOut_ = new double[dimY];

        int myRank;
        MPI_Comm_rank(MPI_COMM_WORLD, &myRank);
        myGridX_ = myRank % numGridsX_;
        myGridY_ = myRank / numGridsX_;

        initGrid(numInitialPerturbations);
    }

    ~Grid() {
        delete[] myGridPrev_;
        delete[] myGridCurr_;
        delete[] myGridNext_;
        delete[] leftEdgeIn_;
        delete[] rightEdgeIn_;
        delete[] leftEdgeOut_;
        delete[] rightEdgeOut_;
    }

    void doIterations(int numIterations) {
        for (int i = 0; i < numIterations; ++i) {
            exchangeEdges();
            doOneIteration();

            double *tmp = myGridPrev_;
            myGridPrev_ = myGridCurr_;
            myGridCurr_ = myGridNext_;
            myGridNext_ = tmp;
        }
    }

  private:
    void initGrid(int numInitialPerturbations) {
        for(int s = 0; s < numInitialPerturbations; s++){    
            // Determine where to place a circle within the interior of the 2-d domain
            int radius = 20+rand() % 30;
            int xcenter = radius + rand() % (dimX_*numGridsX_ - 2*radius);
            int ycenter = radius + rand() % (dimY_*numGridsY_ - 2*radius);

            // Draw the circle
            for(int x = 1; x < dimX_; x++){
                for(int y = 1; y < dimY_; y++){
                    // The coordinate in the global data array (not just in this rank's portion)
                    int globalx = myGridX_*dimX_ + x;
                    int globaly = myGridY_*dimY_ + y;

                    double distanceToCenter = sqrt((globalx-xcenter)*(globalx-xcenter)
                                                 + (globaly-ycenter)*(globaly-ycenter));

                    if (distanceToCenter < radius) {
                        // ranges from 0 to 3pi/2 
                        double rscaled = (distanceToCenter/radius) * 3.0 * 3.14159/2.0;

                        // Range won't exceed -700 to 700
                        double t = 700.0 * cos(rscaled);

                        myGridCurr_[LIN2(x,y,dimX_)] = myGridPrev_[LIN2(x,y,dimX_)] = t;
                    }
                }						
            }
        }
    }

    void doOneIteration() {
        for (int x = 0; x < dimX_; ++x) {
            for (int y = 0; y < dimY_; ++y) {
                myGridNext_[LIN2(x,y,dimX_)] = 
                    ks * (myGridCurr_[LIN2(x+1,y,dimX_)] + myGridCurr_[LIN2(x-1,y,dimX_)]
                        + myGridCurr_[LIN2(x,y+1,dimX_)] + myGridCurr_[LIN2(x,y-1,dimX_)]
                        - myGridCurr_[LIN2(x,y,dimX_)] * 4)
                  - myGridPrev_[LIN2(x,y,dimX_)] + 2*myGridCurr_[LIN2(x,y,dimX_)];
            }
        }
    }

    void exchangeEdges() {
        MPI_Request reqs[4];
        MPI_Status stats[4];

        int topRank = CYCLIN(myGridX_, myGridY_+1, numGridsX_, numGridsY_);
        int leftRank = CYCLIN(myGridX_-1, myGridY_, numGridsX_, numGridsY_);
        int rightRank = CYCLIN(myGridX_+1, myGridY_, numGridsX_, numGridsY_);
        int bottomRank = CYCLIN(myGridX_, myGridY_-1, numGridsX_, numGridsY_);


        /*************** Recv ***************/
        // Top.
        MPI_Irecv(&(myGridCurr_[LIN(1,0,dimX_+2)]), dimX_, MPI_DOUBLE, topRank, 0, MPI_COMM_WORLD, &reqs[0]);

        // Bottom.
        MPI_Irecv(&(myGridCurr_[LIN(1,dimY_+1,dimX_+2)]), dimX_, MPI_DOUBLE,
                  bottomRank, 0, MPI_COMM_WORLD, &reqs[1]);

        // Left.
        MPI_Irecv(leftEdgeIn_, dimY_, MPI_DOUBLE, leftRank, 0, MPI_COMM_WORLD, &reqs[2]);

        // Right.
        MPI_Irecv(rightEdgeIn_, dimY_, MPI_DOUBLE, rightRank, 0, MPI_COMM_WORLD, &reqs[3]);


        /*************** Send ***************/
        // Top.
        MPI_Send(&(myGridCurr_[LIN2(0,0,dimX_)]), dimX_, MPI_DOUBLE, topRank, 0, MPI_COMM_WORLD);

        // Bottom.
        MPI_Send(&(myGridCurr_[LIN2(0,dimY_-1,dimX_)]), dimX_, MPI_DOUBLE, bottomRank, 0, MPI_COMM_WORLD);

        // Left.
        // A bit more annoying, since the data is not stored contiguously.
        for (int i = 0; i < dimY_; ++i) leftEdgeOut_[i] = myGridCurr_[LIN2(0, i, dimX_)];
        MPI_Send(leftEdgeOut_, dimY_, MPI_DOUBLE, leftRank, 0, MPI_COMM_WORLD);

        // Right.
        for (int i = 0; i < dimY_; ++i) rightEdgeOut_[i] = myGridCurr_[LIN2(dimX_-1, i, dimX_)];
        MPI_Send(rightEdgeOut_, dimY_, MPI_DOUBLE, rightRank, 0, MPI_COMM_WORLD);


        /*************** Update ***************/
        // Wait for all messages to be sent/received.
        MPI_Waitall(4, reqs, stats);

        for (int y = 0; y < dimY_; ++y) {
            myGridCurr_[LIN(0, y+1, dimX_+2)] = leftEdgeIn_[y];
            myGridCurr_[LIN(dimX_+1, y+1, dimX_+2)] = rightEdgeIn_[y];
        }
    }
};

int main(int argc, char *argv[])
{
    // Initialize MPI stuff.
    int nprocs, myRank;
    MPI_Init(&argc, &argv);
    MPI_Comm_size(MPI_COMM_WORLD, &nprocs);
    MPI_Comm_rank(MPI_COMM_WORLD, &myRank);

    // Seed the RNG.
    int seed;
    if (myRank == 0) seed = time(NULL);
    MPI_Bcast(&seed, 1, MPI_INT, 0, MPI_COMM_WORLD);
    srand(seed);

    // Complain if we're missing arguments.
    assert(argc >= 6);

    // Get dimensionality information from the command line.
    // FIXME?: Number of ranks cannot be prime.
    int myDimX = atoi(argv[1]), myDimY = atoi(argv[2]);
    int numGridsX = atoi(argv[3]), numGridsY = atoi(argv[4]);
    int numInitPerturbations = atoi(argv[5]);

    // Make sure those dimensions are right before moving on.
    assert(numGridsX * numGridsY == nprocs);

    Grid grid(myDimX, myDimY, numGridsX, numGridsY, numInitPerturbations);
    grid.doIterations(100);

    MPI_Finalize();
    return 0;
}