/******************************************************************************/ /** @file "MPI/C/Primer 5 - wave/wave_mpi.c" @brief Paralelno resenje 1D Wave Equation @author Milos Gligoric @author Andrija Bosnjakovic */ /******************************************************************************/ #include #include #include /******************************************************************************/ /** @def MIN_NUM_OF_PROCS @brief Minimalni broj procesa u sistemu. */ #define MIN_NUM_OF_PROCS 2 /** @def MAX_NUM_OF_PROCS @brief Maksimalni broj procesa u sistemu. */ #define MAX_NUM_OF_PROCS 5 /** @def MIN_NUM_OF_STEPS @brief Minimalni broj koraka za izracunavanje amplitude. */ #define MIN_NUM_OF_STEPS 1 /** @def MAX_NUM_OF_STEPS @brief Maksimalni broj koraka za izracunavanje amplitude. */ #define MAX_NUM_OF_STEPS 100 /** @def ErrorNumOfProcs @brief Greska u slucaju nekorektnog broja procesa. */ #define ErrorNumOfProcs 2 /** @def ErrorNumOfPoints @brief Greska u slucaju nekorektnog broja tacaka(amplituda) */ #define ErrorNumOfPoints 3 /** @def ErrorNumOfSteps @brief Greska u slucaju nekorektnog broja koraka za izracunavanje amplitude. */ #define ErrorNumOfSteps 4 /** @def ErrorAlloc @brief Greska u slucaju neuspele alokacije. */ #define ErrorAlloc 5 /******************************************************************************/ /** @brief Kreiranje niza zadate duzine @param pLength Duzina niza koji se kreira @return Novokreirani niz */ double* createArray(int pLength) { double *dArray = (double*) calloc(pLength, sizeof(double)); return dArray; } /** @brief Kreiranje niza zadate duzine u slucaju greske u alokaciji poziva Abort @param pLength Duzina niza koji se kreira @return Novokreirani niz */ double* createArrayIfErrorAbort(int pLength) { double *dArray = createArray(pLength); if (dArray == NULL) { MPI_Abort(MPI_COMM_WORLD, ErrorAlloc); } return dArray; } /** @brief Oslobadja memoriju koji je zauzimao dati niz @param pArray Niz koji treba osloboditi */ void destroyArray(double *pArray) { free(pArray); } /** @brief Ucitava zadati broj elemenata i smesta u dati niz @param pArray Prostor u koji treba ucitati elemente @param pLength Broj elemenata koje treba ucitati */ void readArray(double *pArray, int pLength) { int i; for (i = 0; i < pLength; i++) { scanf("%lf", &pArray[i]); } } /** @brief Ispisivanje datog niza uz zadato formatiranje za svaki element @param pArray Niz ciji ce elementi biti ispisani @param pLength Duzina niza koji se ispisuje @param pFormat Format koji se koristi pri ispisu elementa niza */ void printArray(double *pArray, int pLength, const char *pFormat) { int i; for (i = 0; i < pLength; i++) { printf(pFormat, pArray[i]); } } /** @brief Kopira niz sa izvorisne na odredisnu adresu i to zadati broj elemenata @param pSrc Pocetna adresa izvorisnog niza @param pDst Pocetna adresa odredisnog niza @param pLength Broj elemenata koji se kopira */ void copyArray(double *pSrc, double *pDst, int pLength) { int i; for (i = 0; i < pLength; i++) { pDst[i] = pSrc[i]; } } /******************************************************************************/ /** @brief Proverava ispravnost broja procesa @param pNumOfProcs Broj procesa @return 0 ukoliko je broj korektan i !=0 ukoliko je broj nekorektan */ int isWrongNumOfProcs(int pNumOfProcs) { return (pNumOfProcs < MIN_NUM_OF_PROCS || pNumOfProcs > MAX_NUM_OF_PROCS); } /** @brief Ukoliko je broj procesa nekorektan prekida izvrsavanje programa @param pNumOfProcs Broj procesa */ void ifWrongNumOfProcsAbort(int pNumOfProcs) { if (isWrongNumOfProcs(pNumOfProcs)) { MPI_Abort(MPI_COMM_WORLD, ErrorNumOfProcs); } } /** @brief Provera broja tacaka @param pNumOfPoints Broj tacaka @param pSize Broj procesa @return 0 ukoliko je zadati broj korektan i !=0 ukoliko je broj nekorektan */ int isWrongNumOfPoints(int pNumOfPoints, int pSize) { return (pNumOfPoints < pSize || pNumOfPoints%pSize != 0); } /** @brief Ukoliko je broj nekorektan prekida izvrsavanje programa @param pNumOfPoints Broj tacaka @param pSize Broj procesa */ void ifWrongNumOfPointsAbort(int pNumOfPoints, int pSize) { if (isWrongNumOfPoints(pNumOfPoints, pSize)) { MPI_Abort(MPI_COMM_WORLD, ErrorNumOfPoints); } } /** @brief Provera broja koraka za izracunavanje krajnjih amplituda @param pNumOfSteps Broj koraka koje treba izvrsiti @return 0 ukoliko je zadati broj korektan i !=0 ukoliko je broj nekorektan */ int isWrongNumOfSteps(int pNumOfSteps) { return (pNumOfSteps < MIN_NUM_OF_STEPS || pNumOfSteps > MAX_NUM_OF_STEPS); } /** @brief Ukoliko je broj nekorektan prekida izvrsavanje programa @param pNumOfSteps Broj koraka koje treba izvrsiti */ void ifWrongNumOfStepsAbort(int pNumOfSteps) { if (isWrongNumOfSteps(pNumOfSteps)) { MPI_Abort(MPI_COMM_WORLD, ErrorNumOfSteps); } } /** @brief Izracunava vrednost pow(tau, 2) @return pow(tau, 2) */ double sqtau() { double dtime = 0.3; double c = 1; double dx = 1; double tau = (dtime * c / dx); return tau * tau; } /** @brief Izracunava rang levog i desnog suseda @param pRank Rang procesa ciji se susedi odredjuju @param pSize Broj procesa @param pLeft Rang levog procesa @param pRight Rang desnog procesa */ void calculateLeftRightNeighbour(int pRank, int pSize, int *pLeft, int *pRight) { if (pRank == 0) (*pLeft) = pSize - 1; else (*pLeft) = pRank - 1; if (pRank == pSize - 1) (*pRight) = 0; else (*pRight) = pRank + 1; } /******************************************************************************/ /** @brief Ispis zaglavlja @param pRank Rang procesa @param pSize Broj procesa u svetu @param pLeftRank Rang levog suseda @param pRightRank Rang desnog suseda @param pArray Niz podataka koje je proces dobio na obradu @param pLength Duzina niza koji proces obradjuje */ void printHeader(int pRank, int pSize, int pLeftRank, int pRightRank, double *pArray, int pLength) { printf("---------------------------------------------------\n"); printf("#%d od %d procesa\n", pRank, pSize); printf("levi proces #%d\ndesni proces #%d\n", pLeftRank, pRightRank); printf("---------------------------------------------------\n"); printf("iteracija |\n"); printf("---------------------------------------------------\n"); printf("%-#10d|", 0); printArray(pArray, pLength, "%#10g "); printf("\n"); } /** @brief Ispis informacija vezan za jednu iteraciju obrade @param pIteration Trenutna iteracija @param pArray Niz koji opisuje stanje amplituda koje proces obradjuje @param pLength Duzina niza amplituda */ void printIteration(int pIteration, double *pArray, int pLength) { printf("%-#10d|", pIteration); printArray(pArray, pLength, "%#10g "); printf("\n"); } /** @brief Ispis nakon obrade podataka */ void printFooter() { printf("---------------------------------------------------\n"); } /******************************************************************************/ /** @brief Inicijalizacija Master procesa @param pNumOfPoints Broj tacaka (izlazni parametar) @param pPoints Pocetni niz amplituda (izlazni parametar) @param pNumOfSteps Broj koraka koje treba sprovesti (izlazni parametar) @param pSize Broj procesa @param pRank Rang procesa */ void initMaster(int *pNumOfPoints, double **pPoints, int *pNumOfSteps, int pSize, int pRank) { printf("Unesite broj tacaka na x osi: "); scanf("%d", pNumOfPoints); ifWrongNumOfPointsAbort(*pNumOfPoints, pSize); MPI_Bcast(pNumOfPoints, 1, MPI_INT, 0, MPI_COMM_WORLD); /* kreiranje niza amplituda i unos pocetnih vrednosti */ *pPoints = createArrayIfErrorAbort(*pNumOfPoints); printf("Unesite pocetne vrednosti amplituda. Njih %d: ", *pNumOfPoints); readArray(*pPoints, *pNumOfPoints); printf("Unesite broj koraka za izracunavanje amplituda: "); scanf("%d", pNumOfSteps); ifWrongNumOfStepsAbort(*pNumOfSteps); MPI_Bcast(pNumOfSteps, 1, MPI_INT, 0, MPI_COMM_WORLD); } /** @brief Zavrsne radnje Master procesa @param pNumOfPoints Broj tacaka @param pPoints Pocetni niz amplituda @param pNumOfSteps Broj koraka koje treba sprovesti @param pSize Broj procesa @param pRank Rang procesa */ void finishMaster(int pNumOfPoints, double *pPoints, int pNumOfSteps, int pSize, int pRank) { printf("*** Konacni izgled amplituda ***\n"); printArray(pPoints, pNumOfPoints, "%#10g "); destroyArray(pPoints); } /** @brief Inicijalizacija ostalih procesa @param pNumOfPoints Broj tacaka (izlazni parametar) @param pPoints Pocetni niz amplituda (izlazni parametar) @param pNumOfSteps Broj koraka koje treba sprovesti (izlazni parametar) @param pSize Broj procesa @param pRank Rang procesa */ void initWorker(int *pNumOfPoints, double **pPoints, int *pNumOfSteps, int pSize, int pRank) { MPI_Bcast(pNumOfPoints, 1, MPI_INT, 0, MPI_COMM_WORLD); MPI_Bcast(pNumOfSteps, 1, MPI_INT, 0, MPI_COMM_WORLD); *pPoints = createArrayIfErrorAbort(*pNumOfPoints); } /** @brief Zavrsne radnje Procesa @param pNumOfPoints Broj tacaka @param pPoints Pocetni niz amplituda @param pNumOfSteps Broj koraka koje treba sprovesti @param pSize Broj procesa @param pRank Rang procesa */ void finishWorker(int pNumOfPoints, double *pPoints, int pNumOfSteps, int pSize, int pRank) { destroyArray(pPoints); } /** @brief Posao koji procesi trebaju da obave @param pNumOfPoints Broj tacaka @param pPoints Pocetni niz amplituda @param pNumOfSteps Broj koraka koje treba sprovesti @param pSize Broj procesa @param pRank Rang procesa */ void work(int pNumOfPoints, double *pPoints, int pNumOfSteps, int pSize, int pRank) { int tag = 123; double *myPoints, *oldValues, *newValues; int myNumOfPoints; int leftRank, rightRank; int i, j; /* odredjivanje suseda */ calculateLeftRightNeighbour(pRank, pSize, &leftRank, &rightRank); myNumOfPoints = pNumOfPoints/pSize; myPoints = createArrayIfErrorAbort(myNumOfPoints); /* slanje/prijem podataka za obradu */ MPI_Scatter(pPoints, myNumOfPoints, MPI_DOUBLE, myPoints, myNumOfPoints, MPI_DOUBLE, 0, MPI_COMM_WORLD); oldValues = createArrayIfErrorAbort(myNumOfPoints); copyArray(myPoints, oldValues, myNumOfPoints); newValues = createArrayIfErrorAbort(myNumOfPoints); printHeader(pRank, pSize, leftRank, rightRank, myPoints, myNumOfPoints); /* glavna obrada */ for (i = 0; i < pNumOfSteps; i++) { double leftValue, rightValue; double *temp; MPI_Status status; MPI_Request request; MPI_Isend(&myPoints[myNumOfPoints-1], 1, MPI_DOUBLE, rightRank, tag, MPI_COMM_WORLD, &request); MPI_Isend(&myPoints[0], 1, MPI_DOUBLE, leftRank, tag, MPI_COMM_WORLD, &request); MPI_Recv(&leftValue, 1, MPI_DOUBLE, leftRank, tag, MPI_COMM_WORLD, &status); MPI_Recv(&rightValue, 1, MPI_DOUBLE, rightRank, tag, MPI_COMM_WORLD, &status); for (j = 0; j < myNumOfPoints; j++) { if (j != 0) leftValue = myPoints[j-1]; if (j != myNumOfPoints-1) rightValue = myPoints[j+1]; newValues[j] = (2 * myPoints[j]) - oldValues[j] + (sqtau() * (leftValue - (2 * myPoints[j]) + rightValue)); } temp = newValues; newValues = oldValues; oldValues = myPoints; myPoints = temp; printIteration(i+1, myPoints, myNumOfPoints); } printFooter(); // grupisanje svih amplituda u jedan niz MPI_Gather(myPoints, myNumOfPoints, MPI_DOUBLE, pPoints, myNumOfPoints, MPI_DOUBLE, 0, MPI_COMM_WORLD); destroyArray(myPoints); destroyArray(oldValues); destroyArray(newValues); } /******************************************************************************/ /** @brief Glavni program @param argc Broj argumenata komandne linije @param argv Argumenti komandne linije @return Vraca @b 0 ako je sve uspesno izvrseno */ int main(int argc, char *argv[]) { // rang i velicina MPI sveta int rank, size; /* ukupni broj amplituda i broj amplituda po procesu */ int numOfPoints; /* amplitude zice */ double *points; /* broj koraka za racunanje amplituda */ int numOfSteps; /* inicijalizujemo MPI */ MPI_Init(&argc, &argv); MPI_Comm_rank(MPI_COMM_WORLD, &rank); MPI_Comm_size(MPI_COMM_WORLD, &size); ifWrongNumOfProcsAbort(size); if (rank == 0) { initMaster(&numOfPoints, &points, &numOfSteps, size, rank); } else { initWorker(&numOfPoints, &points, &numOfSteps, size, rank); } work(numOfPoints, points, numOfSteps, size, rank); if (rank == 0) { finishMaster(numOfPoints, points, numOfSteps, size, rank); } else { finishWorker(numOfPoints, points, numOfSteps, size, rank); } /* Zavrsavamo MPI */ MPI_Finalize(); return 0; } /******************************************************************************/