例題(2次元衝撃波管問題)並列プログラムの解説

 

並列プログラムは参考文献[2]の並列プログラミング虎の巻MPI版を参考にしています。

プログラムの解説は、参考文献[2]の該当ページを示します。それ以外のことについてはコメントを記します。

例題のプログラムは動的負荷分散型で設計しています。使用したパフォーマンスファイルはこちらに示しています

 

 

2次元衝撃波管問題並列プログラム


 

 


#include <stdio.h>

#include <math.h>

#include <mpi.h>

#include <time.h>

 

/* Definition of constants */

#define NMAX 1000             シミュレーション回数

#define JMAX 30    y方向のメッシュ数

#define IMAX 20    x方向のメッシュ数

 

/* Definition of Parallel constants */

#define Max_CPU 4              プロセッサ数

#define X_div 4      x方向の分割数

#define Y_div 1       y方向の分割数

 

#define X_MIN (1)

#define Y_MIN (1)

#define X_MAX (IMAX-1)

#define Y_MAX (JMAX-1)

 

↓動的負荷分散型でのみ使用

#define MEASURE_LOOP1 (100000)    プロセッサ処理能力測定ループ繰り返し数

#define MEASURE_LOOP2 (100000)  プロセッサ処理能力測定ループ繰り返し数

#define BASIS (0.200000)     処理時間の変動基準(20%)

#define JUDGEMENT (500)  処理能力を判断する間隔500step

 

/* prototypes definition of all functions in this program */

int main();

void SlvDelta_U();  ΔUのソルバーを含む関数

void SlvDelta_V();  ΔV

void SlvP();            P

void SlvRHO();       RHO(ρ)

void SlvQx();         x方向の人口粘性Qx

void SlvQy();         y方向の人口粘性Qy

void SlvQ();           Q

void SlvT();           T

void SlvU();           U

void SlvV();           V

void MPIinit();        MPI起動

void Set_Parameter();            並列計算をするためのパラメータ設定関数

void Set_Gather();  シミュレーションデータをプロセス0へ収集するためのパラメータ設定関数

void init();              変数を初期化する関数

void M_AVS_out(); MicroAVS用のデータファイル出力関数

void M_AVS_Fld_out();         MicroAVS用のフィールドファイルヘッダー出力

void M_AVS_Fld_out2();       MicroAVS用のフィールドファイル出力

void M_POST_out();             Micro/Post用のデータファイル出力

void Data_Shift();   境界データ通信関数

void Data_Shift2(); 境界より一つ深い位置のデータ通信関数

void Data_Gather();               シミュレーションデータをプロセス0へ収集するための関数

void output();         計算データを標準出力する関数

void subst();           時刻n+1nへの変数渡しを行う関数

void MPIfinish();    MPI終了

void Set_Range();   分割領域のメッシュ幅を設定する関数

void Para_Range(); 領域を均等分割するする関数

void Para_Type_Block2d();    派生データタイプの設定

void Para_Range_Static();      領域を処理能力値によって分割する関数

void Perform_Total();            能力値を総和する関数

void Mesh_Resize();              メッシュ幅をリサイズする関数

void Judge();          プロセッサの処理能力の変化を判断する関数

void CPU_Measure();            プロセッサの処理能力を測定する関数

void Time_AllGather();          実行時間のデータをすべてのプロセスへ送信する関数

void Data_AllGather();           シミュレーションデータをすべてのプロセスへ送信する関数

void Time_Ranking();            処理能力の高い順にプロセスをランク付けする関数

void Reset_Parameter();         並列計算のためのパラメータを再設定する関数

 

/*  Declare the global variables used in Domain Divided Parallel computation  */

/* MPI type Variables */ MPI変数 

MPI_Status status;

MPI_Request isend1;

MPI_Request isend2;

MPI_Request jsend1;

MPI_Request jsend2;

MPI_Request irecv1;

MPI_Request irecv2;

MPI_Request jrecv1;

MPI_Request jrecv2;

MPI_Request data_send;

MPI_Request data_recv[Max_CPU];

MPI_Request data_all_send[Max_CPU];

MPI_Request data_all_recv[Max_CPU];

MPI_Datatype NEW_DATA_TYPE[Max_CPU];

 

/* Int type Variables */

int myrank;            自分のプロセス番号

int nprocs;             プロセス数(並列数)

int myranki;            x方向のプロセス番号

int istart;  分割領域のx方向の始端メッシュ番号

int iend;   分割領域のx方向の終端メッシュ番号

int istart2;              分割領域のx方向の始端メッシュ番号

int iend2; 分割領域のx方向の終端メッシュ番号

int iup;     自分のプロセスより右に位置するプロセス番号

int idown;              自分のプロセスより左に位置するプロセス番号

int ilen;    担当するx方向のメッシュ長さ

int myrankj; y方向の〃

int jstart;

int jend;

int jstart2;

int jend2;

int jup;

int jdown;

int jlen;

int table[X_div+2][Y_div+2];  プロセス番号のマップ配列

int x_mesh_size[X_div];        x方向の各分割領域のメッシュ数

int y_mesh_size[Y_div];         y方向の各分割領域のメッシュ数

int perform_rank[Max_CPU] = {0, 2, 3, 1 };        処理能力のランクデータ

int start_flg;           フラグ

 

/* Double type Variables */

double total;

double GCLOCK;    CLOCK_PER_SECの逆数

double prev_perform;            500stepの計算時間

double calc_perform[Max_CPU];          各プロセスの計算能力値格納配列

double perform[Max_CPU] = {1.000000, 0.490000, 0.970000, 0.960000 };       パフォーマンスファイルからのデータ、各プロセスの計算能力値

 

/* Time type Variables */

time_t start_time;    測定開始変数

time_t end_time;     測定終了変数

clock_t start_clock;               測定開始クロック変数

clock_t end_clock; 測定終了クロック変数

 

/* time step variable */

int n;

/* Micro AVS output counter variable */

int avs_cnt;            MicroAVSカウント変数

 

/* Declaration of parameters */

double Cv = 1.000000;

double m = 1.000000;

double K = 1.000000;

double c = 1.000000;

 

/* Declaration of independent variables */

double dx = 1.000000;           メッシュ幅

double dy = 1.000000;           メッシュ幅

double dt = 0.100000;            時間刻み幅

 

/* Declaration of dependent variables */

double U[IMAX+1][JMAX+1];/* U n timestep = 0  0 */ 変数Un

double U1[IMAX+1][JMAX+1];/* U n timestep = 1  1 */   変数Un+1

double V[IMAX+1][JMAX+1];/* V n timestep = 0  0 */ 変数Vn

double V1[IMAX+1][JMAX+1];/* V n timestep = 1  1 */   変数Vn+1

double RHO[IMAX+1][JMAX+1];/* RHO n timestep = 0  0 */

double RHO1[IMAX+1][JMAX+1];/* RHO n timestep = 1  1 */

double Qx[IMAX+1][JMAX+1];/* Qx n timestep = 0  0 */

double Qx1[IMAX+1][JMAX+1];/* Qx ~ timestep = 0  1 */

double Qy[IMAX+1][JMAX+1];/* Qy n timestep = 0  0 */

double Qy1[IMAX+1][JMAX+1];/* Qy ~ timestep = 0  1 */

double Delta_U[IMAX+1][JMAX+1];/* Delta_U n timestep = 0  0 */

double Delta_U1[IMAX+1][JMAX+1];/* Delta_U ~ timestep = 0  1 */

double Delta_V[IMAX+1][JMAX+1];/* Delta_V n timestep = 0  0 */

double Delta_V1[IMAX+1][JMAX+1];/* Delta_V ~ timestep = 0  1 */

double P[IMAX+1][JMAX+1];/* P n timestep = 0  0 */

double P1[IMAX+1][JMAX+1];/* P ~ timestep = 0  1 */

double Q[IMAX+1][JMAX+1];/* Q n timestep = 0  0 */

double Q1[IMAX+1][JMAX+1];/* Q ~ timestep = 0  1 */

double T[IMAX+1][JMAX+1];/* T n timestep = 0  0 */

double T1[IMAX+1][JMAX+1];/* T n timestep = 1  1 */

 

 

/* Declaration of FILE Pointer */

FILE *filePnt1;

FILE *filePnt2;

FILE *filePnt3;

 

 

int main(int argc,char **argv)

{

              filePnt1 = fopen("M_AVS.dat", "w"); MicroAVSファイル

              filePnt2 = fopen("M_POST.dat", "w");   Micro/POSTファイル

              filePnt3 = fopen("M_AVS.fld", "w");      MicroAVS.fldファイル

 

              MPIinit(argc,argv); MPIの起動

              Set_Parameter();     領域分割のためのパラメータの設定

              Set_Gather();         データ収集のためのパラメータ設定

              init();       変数の初期化

              M_AVS_Fld_out();  フィールドファイルにヘッダー書き出し

 

              /* step 0 calc. */ ステップ0の計算

              {

                             Data_Shift(U,ilen,jlen); 変数Uの分割領域境界データ通信

                             Data_Shift2(U,ilen,jlen);         変数Uの分割領域境界より一つ深い位置のデータ通信

                             SlvDelta_U();         変数Uの計算

                             Data_Shift(V,ilen,jlen);           変数Vの境界データ通信

                             Data_Shift2(V,ilen,jlen);         変数Vの境界より一つ深い位置のデータ通信

                             SlvDelta_V();

                             SlvP();

                             Data_Shift(RHO,ilen,jlen);

                             SlvRHO();

                             SlvQx();

                             SlvQy();

                             Data_Shift2(RHO,ilen,jlen);

                             SlvQ();

                             Data_Shift(T,ilen,jlen);

                             SlvT();

                             Data_Shift(P1,ilen,jlen);

                             Data_Shift(Q1,ilen,jlen);

                             SlvU();

                             SlvV();

              }

              GCLOCK = 1.000000/(CLOCKS_PER_SEC); GCLOCKの計算

              start_flg = 1;

              start_clock = clcok();

              time(&start_time);  時間測定開始

              for(n=1; n<=NMAX; n++) { シミュレーション計算開始、NMAX回繰り返す

                             /* 1 time step solver n --> n+1 */

                             {

                                           Data_Shift(U,ilen,jlen);

                                           Data_Shift2(U,ilen,jlen);

                                           SlvDelta_U();

                                           Data_Shift(V,ilen,jlen);

                                           Data_Shift2(V,ilen,jlen);

                                           SlvDelta_V();

                                           SlvP();

                                           Data_Shift(RHO,ilen,jlen);

                                           SlvRHO();

                                           SlvQx();

                                           SlvQy();

                                           Data_Shift2(RHO,ilen,jlen);

                                           SlvQ();

                                           Data_Shift(T,ilen,jlen);

                                           SlvT();

                                           Data_Shift(P1,ilen,jlen);

                                           Data_Shift(Q1,ilen,jlen);

                                           SlvU();

                                           SlvV();

                             }

                             if((n%20)==0) { step20ごとにデータをプロセス0へ収集し、ファイル出力する

                                           Data_Gather(Delta_U1); 変数ΔUn+1のデータをプロセス0へ送信する

                                           Data_Gather(Delta_V1);

                                           Data_Gather(P1);

                                           Data_Gather(RHO1);

                                           Data_Gather(Qx1);

                                           Data_Gather(Qy1);

                                           Data_Gather(Q1);

                                           Data_Gather(T1);

                                           Data_Gather(U1);

                                           Data_Gather(V1);

                                           if(myrank==0) { プロセス0の処理

                                                         output();  ファイル出力

                                                         M_AVS_out();        MicroAVSデータファイル出力

                                                         M_AVS_Fld_out2();              MicroAVSフィールドファイル出力

                                                         M_POST_out();     Micro/Postデータファイル出力

                                           }

                             }

                             if((n%JUDGEMENT)==0) { 500step毎に

                                           end_clock = clcok();

                                           time(&end_time);    時間測定終了

                                           Judge();   プロセッサの処理能力の変動を判断する

                                           start_clock = clcok();           

                                           time(&start_time);  時間計測開始

                             }

                             subst();   時間n+1の変数を時間nの変数へ代入する

              }

              MPIfinish();           MPIの終了

              return(0); プログラムの終了

              fclose(filePnt1);      ファイルクローズ

              fclose(filePnt2);      ファイルクローズ

              fclose(filePnt3);      ファイルクローズ

 

}

 

 

void SlvDelta_U()   変数ΔUn → ΔUn+1ソルバー

{

              int i, j;

              for(i=istart; i<=iend2; i++) {

                             for(j=jstart; j<=jend; j++) {

                                           Delta_U1[i][j] =  0.500*(U[i+2][j]+U[i+1][j])

                                                            - 0.500*(U[i+1][j]+U[i][j]);

                             }

              }

}

 

 

void SlvDelta_V()   変数ΔVn → ΔVn+1ソルバー

{

              int i, j;

              for(i=istart; i<=iend; i++) {

                             for(j=jstart; j<=jend2; j++) {

                                           Delta_V1[i][j] =  0.500*(V[i][j+2]+V[i][j+1])

                                                            - 0.500*(V[i][j+1]+V[i][j]);

                             }

              }

}

 

 

void SlvP()             変数Pn → Pn+1ソルバー

{

              int i, j;

              for(i=istart; i<=iend; i++) {

                             for(j=jstart; j<=jend; j++) {

                                           P1[i][j] = K/(m)*RHO[i][j]*T[i][j];

                             }

              }

}

 

 

void SlvRHO()        変数ρn → ρn+1ソルバー

{

              int i, j;

              for(i=istart2; i<=iend; i++) {

                             for(j=jstart2; j<=jend; j++) {

                                           RHO1[i][j] =  1.000/( 1.000/(dt))

                                                        *(-(-RHO[i][j])/(dt)

                                                          +(-(RHO[i][j]*(U[i+1][j]-U[i][j])/(dx)

                                                              +( 0.500*(U[i+1][j]+U[i][j])

                                                                +fabs( 0.500*(U[i+1][j]+U[i][j])))/( 2.000)

                                                               *(RHO[i][j]-RHO[i-1][j])/(dx)

                                                              +( 0.500*(U[i+1][j]+U[i][j])

                                                                -fabs( 0.500*(U[i+1][j]+U[i][j])))/( 2.000)

                                                               *(RHO[i+1][j]-RHO[i][j])/(dx))

                                                            -(RHO[i][j]*(V[i][j+1]-V[i][j])/(dy)

                                                              +( 0.500*(V[i][j+1]+V[i][j])

                                                                +fabs( 0.500*(V[i][j+1]+V[i][j])))/( 2.000)

                                                               *(RHO[i][j]-RHO[i][j-1])/(dy)

                                                              +( 0.500*(V[i][j+1]+V[i][j])

                                                                -fabs( 0.500*(V[i][j+1]+V[i][j])))/( 2.000)

                                                               *(RHO[i][j+1]-RHO[i][j])/(dy))));

                             }

              }

              /* Boundary Conditions */

              {

                             /* boundary condition on LEFT */

                             i = 1;

                             for(j=1; j<=JMAX-1; j++) {

                                           RHO1[i][j] = RHO1[i+1][j];

                             }

                             /* boundary condition on RIGHT */

                             i = IMAX-1;

                             for(j=1; j<=JMAX-1; j++) {

                                           RHO1[i][j] = RHO1[i-1][j];

                             }

                             /* boundary condition on BTTM */

                             j = 1;

                             for(i=1; i<=IMAX-1; i++) {

                                           RHO1[i][j] = RHO1[i][j+1];

                             }

                             /* boundary condition on TOP */

                             j = JMAX-1;

                             for(i=1; i<=IMAX-1; i++) {

                                           RHO1[i][j] = RHO1[i][j-1];

                             }

              }

}

 

 

void SlvQx()          変数Qxn → Qxn+1ソルバー

{

              int i, j;

              for(i=istart; i<=iend; i++) {

                             for(j=jstart; j<=jend; j++) {

                                           Qx1[i][j] = c*c*(fabs(Delta_U1[i][j])-Delta_U1[i][j])/( 2.000)

                                                       *(fabs(Delta_U1[i][j])-Delta_U1[i][j])/( 2.000);

                             }

              }

}

 

 

void SlvQy()          変数Qyn → Qyn+1ソルバー

{

              int i, j;

              for(i=istart; i<=iend; i++) {

                             for(j=jstart; j<=jend; j++) {

                                           Qy1[i][j] = c*c*(fabs(Delta_V1[i][j])-Delta_V1[i][j])/( 2.000)

                                                       *(fabs(Delta_V1[i][j])-Delta_V1[i][j])/( 2.000);

                             }

              }

}

 

 

void SlvQ()            変数Qn → Qn+1ソルバー

{

              int i, j;

              for(i=istart; i<=iend; i++) {

                             for(j=jstart; j<=jend; j++) {

                                           Q1[i][j] = (Qx1[i][j]+Qy1[i][j])*RHO[i][j];

                             }

              }

}

 

 

void SlvT()            変数Tn → Tn+1ソルバー

{

              int i, j;

              for(i=istart2; i<=iend2; i++) {

                             for(j=jstart2; j<=jend2; j++) {

                                           T1[i][j] =  1.000/( 1.000/(dt))

                                                      *(-(-T[i][j])/(dt)

                                                        -( 0.500*(U[i+1][j]+U[i][j])

                                                          +fabs( 0.500*(U[i+1][j]+U[i][j])))/( 2.000)

                                                         *(T[i][j]-T[i-1][j])/(dx)

                                                        -( 0.500*(U[i+1][j]+U[i][j])

                                                          -fabs( 0.500*(U[i+1][j]+U[i][j])))/( 2.000)

                                                         *(T[i+1][j]-T[i][j])/(dx)

                                                        -( 0.500*(V[i][j+1]+V[i][j])

                                                          +fabs( 0.500*(V[i][j+1]+V[i][j])))/( 2.000)

                                                         *(T[i][j]-T[i][j-1])/(dy)

                                                        -( 0.500*(V[i][j+1]+V[i][j])

                                                          -fabs( 0.500*(V[i][j+1]+V[i][j])))/( 2.000)

                                                         *(T[i][j+1]-T[i][j])/(dy)

                                                        +(-(P1[i][j]+Q1[i][j])/(Cv*RHO[i][j]))

                                                         *((U[i+1][j]-U[i][j])/(dx)+(V[i][j+1]-V[i][j])/(dy)));

                             }

              }

              /* Boundary Conditions */

              {

                             /* boundary condition on LEFT */

                             i = 1;

                             for(j=1; j<=JMAX-1; j++) {

                                           T1[i][j] = T1[i+1][j];

                             }

                             /* boundary condition on RIGHT */

                             i = IMAX-1;

                             for(j=1; j<=JMAX-1; j++) {

                                           T1[i][j] = T1[i-1][j];

                             }

                             /* boundary condition on BTTM */

                             j = 1;

                             for(i=1; i<=IMAX-1; i++) {

                                           T1[i][j] = T1[i][j+1];

                             }

                             /* boundary condition on TOP */

                             j = JMAX-1;

                             for(i=1; i<=IMAX-1; i++) {

                                           T1[i][j] = T1[i][j-1];

                             }

              }

}

 

 

void SlvU()            変数Un → Un+1ソルバー

{

              int i, j;

              for(i=istart2; i<=iend2; i++) {

                             for(j=jstart2; j<=jend2; j++) {

                                           U1[i][j] =  1.000/( 1.000/(dt))

                                                      *(-(-U[i][j])/(dt)

                                                        -(U[i][j]+fabs(U[i][j]))/( 2.000)

                                                         *(U[i][j]-U[i-1][j])/(dx)

                                                        -(U[i][j]-fabs(U[i][j]))/( 2.000)

                                                         *(U[i+1][j]-U[i][j])/(dx)

                                                        -( 0.250*(V[i][j+1]+V[i][j]+V[i-1][j+1]+V[i-1][j])

                                                          +fabs( 0.250

                                                                *(V[i][j+1]+V[i][j]+V[i-1][j+1]+V[i-1][j])))/( 2.000)

                                                         *(U[i][j]-U[i][j-1])/(dy)

                                                        -( 0.250*(V[i][j+1]+V[i][j]+V[i-1][j+1]+V[i-1][j])

                                                          -fabs( 0.250

                                                                *(V[i][j+1]+V[i][j]+V[i-1][j+1]+V[i-1][j])))/( 2.000)

                                                         *(U[i][j+1]-U[i][j])/(dy)

                                                        +(- 1.000/( 0.500*(RHO[i][j]+RHO[i-1][j])))

                                                         *((P1[i][j]-P1[i-1][j])/(dx)

                                                           +(Q1[i][j]-Q1[i-1][j])/(dx)));

                             }

              }

              /* Boundary Conditions */

              {

                            /* boundary condition on LEFT */

                             i = 1;

                             for(j=1; j<=JMAX-1; j++) {

                                           U1[i][j] = U1[i+2][j]+U1[i+1][j]+U1[i+1][j];

                             }

                             /* boundary condition on RIGHT */

                             i = IMAX-1;

                             for(j=1; j<=JMAX-1; j++) {

                                           U1[i+1][j] = U1[i][j]+U1[i-1][j]+U1[i][j];

                             }

                             /* boundary condition on BTTM */

                             j = 1;

                             for(i=1; i<=IMAX-1; i++) {

                                           U1[i][j] = U1[i][j+1];

                             }

                             /* boundary condition on TOP */

                             j = JMAX-1;

                             for(i=1; i<=IMAX-1; i++) {

                                           U1[i][j] = U1[i][j-1];

                             }

              }

}

 

 

void SlvV()            変数Vn → Vn+1ソルバー

{

              int i, j;

              for(i=istart2; i<=iend2; i++) {

                             for(j=jstart2; j<=jend2; j++) {

                                           V1[i][j] =  1.000/( 1.000/(dt))

                                                      *(-(-V[i][j])/(dt)

                                                        -( 0.250*(U[i+1][j]+U[i+1][j-1]+U[i][j]+U[i][j-1])

                                                          +fabs( 0.250

                                                                *(U[i+1][j]+U[i+1][j-1]+U[i][j]+U[i][j-1])))/( 2.000)

                                                         *(V[i][j]-V[i-1][j])/(dx)

                                                        -( 0.250*(U[i+1][j]+U[i+1][j-1]+U[i][j]+U[i][j-1])

                                                          -fabs( 0.250

                                                                *(U[i+1][j]+U[i+1][j-1]+U[i][j]+U[i][j-1])))/( 2.000)

                                                         *(V[i+1][j]-V[i][j])/(dx)

                                                        -(V[i][j]+fabs(V[i][j]))/( 2.000)

                                                         *(V[i][j]-V[i][j-1])/(dy)

                                                        -(V[i][j]-fabs(V[i][j]))/( 2.000)

                                                         *(V[i][j+1]-V[i][j])/(dy)

                                                        +(- 1.000/( 0.500*(RHO[i][j]+RHO[i][j-1])))

                                                         *((P1[i][j]-P1[i][j-1])/(dy)

                                                           +(Q1[i][j]-Q1[i][j-1])/(dy)));

                             }

              }

              /* Boundary Conditions */

              {

                             /* boundary condition on LEFT */

                             i = 1;

                             for(j=1; j<=JMAX-1; j++) {

                                           V1[i][j] = V1[i+1][j];

                             }

                             /* boundary condition on RIGHT */

                             i = IMAX-1;

                             for(j=1; j<=JMAX-1; j++) {

                                           V1[i][j] = V1[i-1][j];

                             }

                             /* boundary condition on BTTM */

                             j = 1;

                             for(i=1; i<=IMAX-1; i++) {

                                           V1[i][j] = V1[i][j+2]+V1[i][j+1]+V1[i][j+1];

                             }

                             /* boundary condition on TOP */

                             j = JMAX-1;

                             for(i=1; i<=IMAX-1; i++) {

                                           V1[i][j+1] = V1[i][j]+V1[i][j-1]+V1[i][j];

                             }

              }

}

 

 

void MPIinit(int char_num,char **char_str) MPIの起動

{             虎の巻 ページ3-2参照

              MPI_Init(&char_num,&char_str);

              MPI_Comm_size(MPI_COMM_WORLD,&nprocs);

              MPI_Comm_rank(MPI_COMM_WORLD,&myrank);

              if(nprocs!=X_div*Y_div) {

                             printf("ERROR: It is difference between nprocs(on MPI environment) and divison number(on P-NCAS)\n");

                             MPI_Finalize();

                             exit(-1);

              }

}

 

 

void Set_Parameter()             領域分割のためのパラメータ設定

{             虎の巻 ページ5-6参照

              int i, j, rank, c;

              c = 1;

              for(j=-1; j<=Y_div; j++) {

                             for(i=-1; i<=X_div; i++) {      プロセステーブルの初期化

                                           table[i+c][j+c] = MPI_PROC_NULL;

                             }

              }

              rank = 0;

              for(j=0; j<=Y_div-1; j++) {

                             for(i=0; i<=X_div-1; i++) {    プロセステーブルの設定

                                           table[i+c][j+c] = rank;

                                           if(myrank==rank) {

                                                         myranki = i;

                                                         myrankj = j;

                                           }

                                           rank = rank+1;

                             }

              }

              if(X_div!=1) {

              Para_Range_Static(X_MIN,X_MAX,X_div,myranki,&istart,&iend,0,x_mesh_size);

              } else {    虎の巻 ページ4-26参照

                             Para_Range(X_MIN,X_MAX,X_div,myranki,&istart,&iend);

              }

              if(Y_div!=1) {

                             Para_Range_Static(Y_MIN,Y_MAX,Y_div,myrankj,&jstart,&jend,1,y_mesh_size);

              } else {

                             Para_Range(Y_MIN,Y_MAX,Y_div,myrankj,&jstart,&jend);

              }

              ilen = iend-istart+1;

              jlen = jend-jstart+1;

              Set_Range(istart,iend,myranki,X_div,&istart2,&iend2);

              Set_Range(jstart,jend,myrankj,Y_div,&jstart2,&jend2);

              iup = table[myranki+1+c][myrankj+c];

              idown = table[myranki-1+c][myrankj+c];

              jup = table[myranki+c][myrankj+1+c];

              jdown = table[myranki+c][myrankj-1+c];

}

 

 

void Set_Gather()   プロセス0へのデータ収集のためのパラメータ設定

{

              int rank, irank, jrank;

              rank = 0;

              for(jrank=0; jrank<=Y_div-1; jrank++) {

                             if(Y_div!=1) {

                                           Para_Range_Static(Y_MIN,Y_MAX,Y_div,jrank,&jstart,&jend,1,y_mesh_size);

                             } else {

                                           Para_Range(Y_MIN,Y_MAX,Y_div,jrank,&jstart,&jend);

                             }

                             for(irank=0; irank<=X_div-1; irank++) {

                                           if(X_div!=1) {

                                                         Para_Range_Static(X_MIN,X_MAX,X_div,irank,&istart,&iend,0,x_mesh_size);

                                           } else {

                                                         Para_Range(X_MIN,X_MAX,X_div,irank,&istart,&iend);

                                           }

              虎の巻ページ3-44参照                     Para_Type_Block2d(1,JMAX,1,istart,iend,jstart,jend,MPI_DOUBLE,&NEW_DATA_TYPE[rank]);

                                           rank = rank+1;

                             }

              }

              if(X_div!=1) {

                             Para_Range_Static(X_MIN,X_MAX,X_div,myranki,&istart,&iend,0,x_mesh_size);

              } else {

                             Para_Range(X_MIN,X_MAX,X_div,myranki,&istart,&iend);

              }

              if(Y_div!=1) {

                             Para_Range_Static(Y_MIN,Y_MAX,Y_div,myrankj,&jstart,&jend,1,y_mesh_size);

              } else {

                             Para_Range(Y_MIN,Y_MAX,Y_div,myrankj,&jstart,&jend);

              }

}

 

 

void init() 各変数の初期化

{

              int i, j;

              for(i=istart; i<=iend; i++) {

                             for(j=jstart; j<=jend; j++) {

                                           Delta_U[i][j] = 0.000000;

                                           Delta_V[i][j] = 0.000000;

                                           Qx[i][j] = 0.000000;

                                           Qy[i][j] = 0.000000;

                                           Q[i][j] = 0.000000;

                                           P[i][j] = 0.000000;

                                           RHO[i][j] = 0.000000;

                                           T[i][j] = 0.000000;

                                           U[i][j] = 0.000000;

                                           V[i][j] = 0.000000;

                             }

              }

}

 

 

void M_AVS_out()

{

              /* Output data file for Micro_AVS */

              int j, i;

              fprintf(filePnt1, "%d \n",n);

              for(j=1; j<=JMAX; j++) {

                             for(i=1; i<=IMAX; i++) {

                                           fprintf(filePnt1, "%d %d %e %e %e %e %e %e %e %e %e %e \n",i,j,Delta_U1[i][j],Delta_V1[i][j],Qx1[i][j],Qy1[i][j],Q1[i][j],P1[i][j],RHO1[i][j],T1[i][j],U1[i][j],V1[i][j]);

                             }

              }

}

 

 

void M_AVS_Fld_out()

{             MicroAVSのマニュアルを参照

              /* Output header field file for Micro_AVS */

              int dim, var_no, skip, stride;

              char filetype[100];

              avs_cnt = 0;

              dim = 2;

              var_no = 10;

              skip = 0;

              stride = var_no+dim;

              strcpy(filetype,"ascii");

              fprintf(filePnt3, "# AVS field file\n");

              fprintf(filePnt3, "ndim = %d \n",dim);

              fprintf(filePnt3, "dim1 = %d \n",IMAX);

              fprintf(filePnt3, "dim2 = %d \n",JMAX);

              fprintf(filePnt3, "nspace = %d \n",dim);

              fprintf(filePnt3, "veclen = %d \n",var_no);

              fprintf(filePnt3, "data = float \n");

              fprintf(filePnt3, "field = uniform \n");

}

 

 

void M_AVS_Fld_out2()

{

              /* Output calculation step field file for Micro_AVS */

              int dim, var_no, skip, offset, stride, cnt, close;

              char filetype[100];

              dim = 2;

              var_no = 10;

              skip = (IMAX*JMAX+1)*avs_cnt;

              avs_cnt = avs_cnt+1;

              close = 1;

              stride = var_no+dim;

              strcpy(filetype,"ascii");

              fprintf(filePnt3, "time file =.\\M_AVS.dat filetype = %s skip = %d close = %d\n",filetype,skip,close);

              skip = skip+1;

              for(cnt=1; cnt<=var_no; cnt++) {

                             offset = cnt+dim-1;

                             fprintf(filePnt3, "variable %d file = .\\M_AVS.dat filetype = %s skip = %d offset = %d stride = %d close = %d \n",cnt,filetype,skip,offset,stride,close);

              }

              fprintf(filePnt3, "EOT \n");

}

 

 

void M_POST_out()

{

              /* Output data file for Micro_POST */

              int i, j;

              for(i=1; i<=IMAX; i++) {

                             for(j=1; j<=JMAX; j++) {

                                           fprintf(filePnt2, "%d %d %e %e %e %e %e %e %e %e %e %e \n",i,j,Delta_U1[i][j],Delta_V1[i][j],Qx1[i][j],Qy1[i][j],Q1[i][j],P1[i][j],RHO1[i][j],T1[i][j],U1[i][j],V1[i][j]);

                             }

              }

}

 

 

void Data_Shift(double Send_Var[IMAX+1][JMAX+1],int i_length,int j_length)  境界データ通信

{             虎の巻ページ5-7参照

              int i, j;

              double Send_Left[JMAX+1], Send_Right[JMAX+1];

              double Send_Top[IMAX+1], Send_Bttm[IMAX+1];

              double Recv_Left[JMAX+1], Recv_Right[JMAX+1];

              double Recv_Top[IMAX+1], Recv_Bttm[IMAX+1];

              境界データのパック

              if(myranki!=X_div-1) { 分割領域のRIGHT境界のパック

                             for(j=jstart; j<=jend; j++) {

                                           Send_Right[j] = Send_Var[iend][j];

                             }

              }

              if(myranki!=0) {     分割領域のLEFT境界のパック

                             for(j=jstart; j<=jend; j++) {

                                           Send_Left[j] = Send_Var[istart][j];

                             }

              }

              if(myrankj!=Y_div-1) {          分割領域のTOP境界のパック

                             for(i=istart; i<=iend; i++) {

                                           Send_Top[i] = Send_Var[i][jend];

                             }

              }

              if(myrankj!=0) {

                             for(i=istart; i<=iend; i++) {    分割領域のBOTTOM境界のパック

                                           Send_Bttm[i] = Send_Var[i][jstart];

                             }

              }

              MPI_Isend(&Send_Right[jstart],j_length,MPI_DOUBLE,iup,1,MPI_COMM_WORLD,&isend1);      RIGHT境界のデータをiup(右隣の)プロセスへ送信

              MPI_Isend(&Send_Left[jstart],j_length,MPI_DOUBLE,idown,1,MPI_COMM_WORLD,&isend2);   LEFT境界のデータをidown(左隣の)プロセスへ送信

              MPI_Isend(&Send_Top[istart],i_length,MPI_DOUBLE,jup,1,MPI_COMM_WORLD,&jsend1);        TOP境界のデータをjup(上隣の)プロセスへ送信

              MPI_Isend(&Send_Bttm[istart],i_length,MPI_DOUBLE,jdown,1,MPI_COMM_WORLD,&jsend2);  BOTTOM境界のデータをjdown(下隣の)プロセスへ送信

              MPI_Irecv(&Recv_Left[jstart],j_length,MPI_DOUBLE,idown,1,MPI_COMM_WORLD,&irecv1);   idown(左隣の)プロセスから送信されてきたデータをRecv_Left配列に受信する

              MPI_Irecv(&Recv_Right[jstart],j_length,MPI_DOUBLE,iup,1,MPI_COMM_WORLD,&irecv2);      iup(右隣の)プロセスから送信されてきたデータをRecv_Right配列に受信する

              MPI_Irecv(&Recv_Bttm[istart],i_length,MPI_DOUBLE,jdown,1,MPI_COMM_WORLD,&jrecv1); jdown(下隣の)プロセスから送信されてきたデータをRecv_Bttm配列に受信する

              MPI_Irecv(&Recv_Top[istart],i_length,MPI_DOUBLE,jup,1,MPI_COMM_WORLD,&jrecv2);        jup(上隣の)プロセスから送信されてきたデータをRecv_Top配列に受信する

              MPI_Wait(&isend1,&status); 同期をとる

              MPI_Wait(&isend2,&status);

              MPI_Wait(&jsend1,&status);

              MPI_Wait(&jsend2,&status);

              MPI_Wait(&irecv1,&status);

              MPI_Wait(&irecv2,&status);

              MPI_Wait(&jrecv1,&status);

              MPI_Wait(&jrecv2,&status);

              境界データのアンパック

              if(myranki!=0) {

                             for(j=jstart; j<=jend; j++) {

                                           Send_Var[istart-1][j] = Recv_Left[j];

                             }

              }

              if(myranki!=X_div-1) {

                             for(j=jstart; j<=jend; j++) {

                                           Send_Var[iend+1][j] = Recv_Right[j];

                             }

              }

              if(myrankj!=0) {

                             for(i=istart; i<=iend; i++) {

                                           Send_Var[i][jstart-1] = Recv_Bttm[i];

                             }

              }

              if(myrankj!=Y_div-1) {

                             for(i=istart; i<=iend; i++) {

                                           Send_Var[i][jend+1] = Recv_Top[i];

                             }

              }

}

 

 

void Data_Shift2(double Send_Var[IMAX+1][JMAX+1],int i_length,int j_length) メッシュ一つ分深い位置のデータ通信

{

              int i, j;

              double Send_Left[JMAX+1], Send_Right[JMAX+1];

              double Send_Top[IMAX+1], Send_Bttm[IMAX+1];

              double Recv_Left[JMAX+1], Recv_Right[JMAX+1];

              double Recv_Top[IMAX+1], Recv_Bttm[IMAX+1];

              if(myranki!=X_div-1) {

                             for(j=jstart; j<=jend; j++) {

                                           Send_Right[j] = Send_Var[iend-1][j];

                             }

              }

              if(myranki!=0) {

                             for(j=jstart; j<=jend; j++) {

                                           Send_Left[j] = Send_Var[istart+1][j];

                             }

              }

              if(myrankj!=Y_div-1) {

                             for(i=istart; i<=iend; i++) {

                                           Send_Top[i] = Send_Var[i][jend-1];

                             }

              }

              if(myrankj!=0) {

                             for(i=istart; i<=iend; i++) {

                                           Send_Bttm[i] = Send_Var[i][jstart+1];

                             }

              }

              MPI_Isend(&Send_Right[jstart],j_length,MPI_DOUBLE,iup,1,MPI_COMM_WORLD,&isend1);

              MPI_Isend(&Send_Left[jstart],j_length,MPI_DOUBLE,idown,1,MPI_COMM_WORLD,&isend2);

              MPI_Isend(&Send_Top[istart],i_length,MPI_DOUBLE,jup,1,MPI_COMM_WORLD,&jsend1);

              MPI_Isend(&Send_Bttm[istart],i_length,MPI_DOUBLE,jdown,1,MPI_COMM_WORLD,&jsend2);

              MPI_Irecv(&Recv_Left[jstart],j_length,MPI_DOUBLE,idown,1,MPI_COMM_WORLD,&irecv1);

              MPI_Irecv(&Recv_Right[jstart],j_length,MPI_DOUBLE,iup,1,MPI_COMM_WORLD,&irecv2);

              MPI_Irecv(&Recv_Bttm[istart],i_length,MPI_DOUBLE,jdown,1,MPI_COMM_WORLD,&jrecv1);

              MPI_Irecv(&Recv_Top[istart],i_length,MPI_DOUBLE,jup,1,MPI_COMM_WORLD,&jrecv2);

              MPI_Wait(&isend1,&status);

              MPI_Wait(&isend2,&status);

              MPI_Wait(&jsend1,&status);

              MPI_Wait(&jsend2,&status);

              MPI_Wait(&irecv1,&status);

              MPI_Wait(&irecv2,&status);

              MPI_Wait(&jrecv1,&status);

              MPI_Wait(&jrecv2,&status);

              if(myranki!=0) {

                             for(j=jstart; j<=jend; j++) {

                                           Send_Var[istart-2][j] = Recv_Left[j];

                             }

              }

              if(myranki!=X_div-1) {

                             for(j=jstart; j<=jend; j++) {

                                           Send_Var[iend+2][j] = Recv_Right[j];

                             }

              }

              if(myrankj!=0) {

                             for(i=istart; i<=iend; i++) {

                                           Send_Var[i][jstart-2] = Recv_Bttm[i];

                             }

              }

              if(myrankj!=Y_div-1) {

                             for(i=istart; i<=iend; i++) {

                                           Send_Var[i][jend+2] = Recv_Top[i];

                             }

              }

}

 

 

void Data_Gather(double Gather_Var[IMAX+1][JMAX+1])  プロセス0へのシミュレーションデータ送信関数

{             虎の巻ページ4-48参照

              int rank;

              if(myrank==0) {

                             for(rank=1; rank<=nprocs-1; rank++) {

                                           MPI_Irecv(Gather_Var,1,NEW_DATA_TYPE[rank],rank,1,MPI_COMM_WORLD,&data_recv[rank]);

                             }

                             for(rank=1; rank<=nprocs-1; rank++) {

                                           MPI_Wait(&data_recv[rank],&status);

                             }

              } else {

                             MPI_Isend(Gather_Var,1,NEW_DATA_TYPE[myrank],0,1,MPI_COMM_WORLD,&data_send);

                             MPI_Wait(&data_send,&status);

              }

}

 

 

void output()

{

              int i, j;

              for(i=1; i<=IMAX; i++) {

                             for(j=1; j<=JMAX; j++) {

                                           printf("%e %e %e %e %e %e %e %e %e %e \n",Delta_U[i][j],Delta_V[i][j],Qx[i][j],Qy[i][j],Q[i][j],P[i][j],RHO[i][j],T[i][j],U[i][j],V[i][j]);

                             }

              }

}

 

 

void subst()            次ステップ計算のため時間項n+1の変数を時間nの変数に入れ替える

{

              int i, j;

              for(i=istart; i<=iend; i++) {

                             for(j=jstart; j<=jend; j++) {

                                           Delta_U[i][j] = Delta_U1[i][j];

                                           Delta_V[i][j] = Delta_V1[i][j];

                                           Qx[i][j] = Qx1[i][j];

                                           Qy[i][j] = Qy1[i][j];

                                           Q[i][j] = Q1[i][j];

                                           P[i][j] = P1[i][j];

                                           RHO[i][j] = RHO1[i][j];

                                           T[i][j] = T1[i][j];

                                           U[i][j] = U1[i][j];

                                           V[i][j] = V1[i][j];

                             }

              }

}

 

 

void MPIfinish()     MPIの終了

{

              MPI_Finalize();

}

 

 

void Set_Range(int start,int end,int axis_rank,int axis_div,int *p_start2,int *p_end2)           計算始端・終端メッシュ番号の設定

{

              *p_start2 = start;

              *p_end2 = end;

              if(axis_rank==0) {

                             *p_start2 = start+1;

              }

              if(axis_rank==axis_div-1) {

                             *p_end2 = end-1;

              }

}

 

 

void Para_Range(int min_mesh,int max_mesh,int proc_num,int rank_num,int *p_start,int *p_end)    領域の均等分割関数

{             虎の巻ページ4-26参照

              int work1, work2, dummy, min;

              work1 = (max_mesh-min_mesh+1)/(proc_num);

              work2 = (max_mesh-min_mesh+1)%proc_num;

              if(work2>=rank_num) {

                             min = rank_num;

              } else {

                             min = work2;

              }

              *p_start = rank_num*work1+min_mesh+min;

              dummy = *p_start;

              *p_end = dummy+work1-1;

              dummy = *p_end;

              if(work2>rank_num) {

                             *p_end = dummy+1;

              }

}

 

 

void Para_Type_Block2d(int cont_min,int cont_max,int discont_min,int i_start,int i_end,int j_start,int j_end,MPI_Datatype old_datatype,MPI_Datatype *p_new_datatype)     派生データタイプの作成

{             虎の巻ページ3-42,44参照

              int i_length, j_length, stride, Data_Size;

              int block[2], disp[2], type[2];

              MPI_Datatype temp_datatype;

              i_length = i_end-i_start+1;

              j_length = j_end-j_start+1;

              MPI_Type_extent(old_datatype,&Data_Size);

              stride = cont_max-cont_min+2;

              MPI_Type_vector(i_length,j_length,stride,old_datatype,&temp_datatype);

              MPI_Type_commit(&temp_datatype);

              block[0] = 1;

              block[1] = 1;

              disp[0] = 0;

              disp[1] = ((cont_max+1)*(i_start-discont_min+1)+(j_start-cont_min+1))

                        *Data_Size;

              type[0] = MPI_LB;

              type[1] = temp_datatype;

              MPI_Type_struct(2,block,disp,type,p_new_datatype);

              MPI_Type_commit(p_new_datatype);

}

 

 

void Para_Range_Static(int min_mesh,int max_mesh,int proc_num,int rank_num,int *p_start,int *p_end,int axis_sw,int axis_mesh_size[Max_CPU])   プロセッサ能力値データから領域を分割する関数

{

              int i, work1;

              Perform_Total();     プロセッサ能力配列の総和を計算する

              Mesh_Resize(min_mesh,max_mesh,proc_num,axis_sw,axis_mesh_size);           メッシュをリサイズする関数

              work1 = 0;

              for(i=0; i<=rank_num; i++) {

                             work1 = work1+axis_mesh_size[i];

              }

              *p_start = work1-(axis_mesh_size[rank_num]-1);

              *p_start = work1;

}

 

 

void Perform_Total()

{

              /* performance total calc. */

              int i;

              total = 0.000000;

              for(i=0; i<=Max_CPU-1; i++) {

                             total = total+perform[i];

              }

}

 

 

void Mesh_Resize(double min_no,double max_no,int proc_num,int axis_sw,int mesh_size[Max_CPU])            メッシュの担当範囲をリサイズする関数

{

              int i, rank, irank, i_store, jrank, j_store, total_mesh, amari, rank, axis_rank, against_proc_num;

              double axis_perform;

              total_mesh = 0;

              if(axis_sw==0) {

                             against_proc_num = Y_div;

              }

              if(axis_sw==1) {

                             against_proc_num = X_div;

              }

              for(i=0; i<=proc_num-1; i++) {

                             axis_perform = 0.000000;

                             for(axis_rank=1; axis_rank<=against_proc_num; axis_rank++) {

                                           if(axis_sw==0) {

                                                         rank = table[i+1][axis_rank];

                                           }

                                           if(axis_sw==1) {

                                                         rank = table[axis_rank][i+1];

                                           }

                                           axis_perform = axis_perform+perform[rank];

                             }

                             mesh_size[i] = (max_no-min_no+1.000000)*axis_perform/(total);

                             total_mesh = total_mesh+mesh_size[i];

              }

              if(total_mesh!=max_no-min_no+1.000000) {

                             amari = (max_no-min_no+1.000000)-total_mesh;

                             for(i=0; i<=amari-1; i++) {

                                           rank = perform_rank[i];

                                           for(irank=1; irank<=X_div; irank++) {

                                                         for(jrank=1; jrank<=Y_div; jrank++) {

                                                                       if(rank==table[irank][jrank]) {

                                                                                      i_store = irank;

                                                                                      j_store = jrank;

                                                                       }

                                                         }

                                           }

                                           if(axis_sw==0) {

                                                         mesh_size[i_store-1] = mesh_size[i_store-1]+1;

                                           }

                                           if(axis_sw==1) {

                                                         mesh_size[j_store-1] = mesh_size[j_store-1]+1;

                                           }

                             }

              }

}

 

 

void Judge()           プロセッサの処理能力が変動したかどうかを判断する関数

{

              int i, Resize_Check, Resize_On;

              double base, all_clock, all_time;

              all_time = difftime(end_time,start_time);

              all_clock = GCLOCK*(end_clock-start_clock);

              if(start_flg==1) {初めは500stepの時間測定のみ

                             start_flg = 0;

                             prev_perform = all_clock/(all_time);     

                             calc_perform[myrank] = prev_perform;               最初の500stepの処理能力値

              } else {

                             calc_perform[myrank] = all_clock/(all_time);       今回の500stepの処理能力値

              }

              Resize_On = 0; リセット

              Resize_Check = 0;  リセット

              if(fabs((calc_perform[myrank]-prev_perform)/(prev_perform))>=BASIS) {      処理能力値が前回500stepの処理能力値に比べBASIS分の変動があれば

                             Resize_Check = 1;  リサイズON

              } else {

                             Resize_Check = 0;  リサイズOFF

              }

              MPI_Allreduce(&Resize_Check,&Resize_On,1,MPI_INT,MPI_SUM,MPI_COMM_WORLD);         すべてのプロセスのResize_Check変数をチェックする

              if(Resize_On>=1) { リサイズONなら

                             start_flg = 1;

                             CPU_Measure(MEASURE_LOOP1,MEASURE_LOOP2);    プロセッサ処理能力測定ルーチンの開始

                             Time_AllGather(perform);      測定時間データをプロセッサ間通信する

                             /* Since reset domain,gaher variable data from all rank. */

                             {             すべての変数データをすべてのプロセスに送受信する

                                           Data_AllGather(Delta_U1);

                                           Data_AllGather(Delta_V1);

                                           Data_AllGather(P1);

                                           Data_AllGather(RHO1);

                                           Data_AllGather(Qx1);

                                           Data_AllGather(Qy1);

                                           Data_AllGather(Q1);

                                           Data_AllGather(T1);

                                           Data_AllGather(U1);

                                           Data_AllGather(V1);

                             }

                             Time_Ranking();     処理能力データをランク付けする

                             Reset_Parameter(); 領域分割によるパラメータをリセットする

                             Set_Gather();         プロセス0へデータを収集するためのパラメータを設定しなおす

              }

}

 

 

void CPU_Measure() プロセッサ処理能力測定関数

{             浮動小数点計算を繰り返しその時間を測定する

              int i, j, int_data;

              double time_diff, double_data;

              time_t loop_sta_time, loop_end_time;

              clock_t loop_sta_clock, loop_end_clock;

              int_data = 100;

              double_data = 1000000.000000;

              MPI_Barrier(MPI_COMM_WORLD);   プロセッサの同期をとる

              loop_sta_clock = clcok();

              time(&loop_sta_time);           時間測定開始

              for(i=0; i<=Loop_Count1; i++) {          浮動小数点計算

                             for(j=0; j<=Loop_Count2; j++) {

                                           int_data = int_data+int_data;

                                           double_data = double_data+double_data;

                                           int_data = int_data*int_data;

                                           double_data = double_data*double_data;

                                           int_data = int_data/(int_data);

                                           double_data = double_data/(double_data);

                             }

              }

              loop_end_clock = clcok();

              time(&loop_end_time);          計測終了

              time_diff = difftime(loop_end_time,loop_sta_time);

              perform[myrank] = GCLOCK*(loop_end_clock-loop_sta_clock)/(time_diff);     perform[自分のプロセス番号]へ処理能力値を格納する

}

 

 

void Time_AllGather(double Send_Array[IMAX+1])            すべてのプロセスにデータを送信する

{             虎の巻ページ3-14参照

              int rank;

              for(rank=0; rank<=nprocs-1; rank++) {

                             MPI_Bcast(&Send_Array[rank],1,MPI_DOUBLE,rank,MPI_COMM_WORLD);

              }

}

 

 

void Data_AllGather(double AllGather_Var[IMAX+1][JMAX+1])         担当している領域のデータをすべてのプロセスへ送信し、担当以外の領域のデータを受信する

{

              int rank;

              for(rank=0; rank<=nprocs-1; rank++) {

                             MPI_Bcast(AllGather_Var,1,NEW_DATA_TYPE[rank],rank,MPI_COMM_WORLD);

              }

}

 

 

void Time_Ranking()             プロセス処理能力値の高い順にランク付けする関数

{

              int i, j, k, JUMP, speed_rank, equal_number;

              int equal_rank[Max_CPU];

              double perform_base;

              for(i=0; i<=Max_CPU-1; i++) {

                             JUMP = 0;

                             speed_rank = 0;

                             equal_number = 0;

                             perform_base = perform[i];

                             for(j=0; j<=Max_CPU-1; j++) {

                                           if(i==j) {

                                                         continue;

                                           }

                                           if(perform_base<perform[j]) {

                                                         speed_rank = speed_rank+1;

                                           } else {

                                                         if(perform_base==perform[j]) {

                                                                       if(i>j) {

                                                                                      JUMP = 1;

                                                                                      break;

                                                                       }

                                                                       equal_rank[0] = i;

                                                                       equal_number = equal_number+1;

                                                                       equal_rank[equal_number] = j;

                                                         }

                                           }

                             }

                             if(JUMP==0) {

                                           if(equal_number==0) {

                                                         perform_rank[speed_rank] = i;

                                           } else {

                                                         for(k=0; k<=equal_number; k++) {

                                                                       perform_rank[speed_rank+k] = equal_rank[k];

                                                         }

                                           }

                             }

              }

}

 

 

void Reset_Parameter()          領域の再分割のためにパラメータを再計算する

{

              if(X_div!=1) {

                             Para_Range_Static(X_MIN,X_MAX,X_div,myranki,&istart,&iend,0,x_mesh_size);

              } else {

                             Para_Range(X_MIN,X_MAX,X_div,myranki,&istart,&iend);

              }

              if(Y_div!=1) {

                             Para_Range_Static(Y_MIN,Y_MAX,Y_div,myrankj,&jstart,&jend,1,y_mesh_size);

              } else {

                             Para_Range(Y_MIN,Y_MAX,Y_div,myrankj,&jstart,&jend);

              }

              ilen = iend-istart+1;

              jlen = jend-jstart+1;

              Set_Range(istart,iend,myranki,X_div,&istart2,&iend2);

              Set_Range(jstart,jend,myrankj,Y_div,&jstart2,&jend2);

}