OCFD_mpi_dev.cpp 17.4 KB
Newer Older
ccfd's avatar
ccfd committed
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
#include "hip/hip_runtime.h"
/*OpenCFD ver 1.4, CopyRight by Li Xinliang, LNM, Institute of Mechanics, CAS, Beijing, Email: lixl@lnm.imech.ac.cn
MPI Subroutines, such as computational domain partation, MPI message send and recv   
只支持N_MSG_SIZE=0, -2  两种通信方式 
*/

#include "mpi.h"
#include "OCFD_mpi.h"
#include "parameters.h"
#include "parameters_d.h"

#include "cuda_commen.h"
#include "cuda_utility.h"
#include "OCFD_mpi_dev.h"


#ifdef __cplusplus
extern "C"{
#endif

// -------------------------------------------------------------------------------------
// Message send and recv at inner boundary (or 'MPI boundary')
void exchange_boundary_xyz_dev(REAL *hostptr , cudaField * devptr)
{
	exchange_boundary_x_dev(hostptr , devptr, Iperiodic[0]);
	exchange_boundary_y_dev(hostptr , devptr, Iperiodic[1]);
	exchange_boundary_z_dev(hostptr , devptr, Iperiodic[2]);
}
// ----------------------------------------------------------------------------------------
void exchange_boundary_x_dev(REAL *hostptr , cudaField * devptr , int Iperiodic1)
{
	exchange_boundary_x_standard_dev(hostptr , devptr, Iperiodic1);
}
// -----------------------------------------------------------------------------------------------
void exchange_boundary_y_dev(REAL *hostptr, cudaField * devptr , int Iperiodic1)
{
	exchange_boundary_y_standard_dev(hostptr , devptr, Iperiodic1);
}
// -----------------------------------------------------------------------------------------------
void exchange_boundary_z_dev(REAL *hostptr , cudaField * devptr , int Iperiodic1)
{
	exchange_boundary_z_standard_dev(hostptr , devptr, Iperiodic1);
}



void exchange_boundary_x_standard_dev(REAL *hostptr , cudaField * devptr, int Iperiodic1)
{
	memcpy_bound_x(hostptr , devptr->ptr , devptr->pitch , D2H , nx_2lap , ny_2lap , nz_2lap);
	exchange_boundary_x_standard(hostptr , Iperiodic1);
	memcpy_bound_x(hostptr , devptr->ptr , devptr->pitch , H2D , nx_2lap , ny_2lap , nz_2lap);
}
// ------------------------------------------------------
void exchange_boundary_y_standard_dev(REAL *hostptr , cudaField * devptr, int Iperiodic1)
{
	memcpy_bound_y(hostptr , devptr->ptr , devptr->pitch , D2H , nx_2lap , ny_2lap , nz_2lap);
	exchange_boundary_y_standard(hostptr , Iperiodic1);
	memcpy_bound_y(hostptr , devptr->ptr , devptr->pitch , H2D , nx_2lap , ny_2lap , nz_2lap);
}
// ------------------------------------------------------------
void exchange_boundary_z_standard_dev(REAL *hostptr, cudaField * devptr , int Iperiodic1)
{
	memcpy_bound_z(hostptr , devptr->ptr , devptr->pitch , D2H , nx_2lap , ny_2lap , nz_2lap);
	exchange_boundary_z_standard(hostptr , Iperiodic1);
	memcpy_bound_z(hostptr , devptr->ptr , devptr->pitch , H2D , nx_2lap , ny_2lap , nz_2lap);
}


/* ===================================================================================================== */


static char mpi_dev_mem_initialized = 0;
static cudaFieldPack *b_xm , *b_xp;
static cudaFieldPack *b_ym , *b_yp;
static cudaFieldPack *b_zm , *b_zp;

void mpi_dev_buffer_attach(){
    // N in bytes
    new_cudaFieldPack(&b_xm , LAP , ny , nz);
    new_cudaFieldPack(&b_xp , LAP , ny , nz);
    new_cudaFieldPack(&b_ym , nx , LAP , nz);
    new_cudaFieldPack(&b_yp , nx , LAP , nz);
    new_cudaFieldPack(&b_zm , nx , ny , LAP);
    new_cudaFieldPack(&b_zp , nx , ny , LAP);
}
void mpi_dev_buffer_detach(){
    delete_cudaFieldPack(b_xm);
    delete_cudaFieldPack(b_xp);
    delete_cudaFieldPack(b_ym);
    delete_cudaFieldPack(b_yp);
    delete_cudaFieldPack(b_zm);
    delete_cudaFieldPack(b_zp);
}
// unsigned int buffer_align_length = 32; // bytes
// void dev_buffer_malloc(unsigned int N){
//     // N in bytes
//     int n = N/buffer_align_length + 1;
//     hipMalloc(n*buffer_align_length);
// }
// void dev_buffer_free(){

// }
// void new_cudaField_buffer(){
    
// }
// void delete_cudaField_buffer(){

// }

void opencfd_mem_init_mpi_dev(){
    if(mpi_dev_mem_initialized == 0){
        mpi_dev_mem_initialized = 1;
        mpi_dev_buffer_attach();
    }
}
void opencfd_mem_finalize_mpi_dev(){
    if(mpi_dev_mem_initialized == 1){
        mpi_dev_mem_initialized = 0;
        mpi_dev_buffer_detach();
    }
}

void exchange_boundary_xyz_packed_dev(REAL *hostptr , cudaField * devptr)
{
	exchange_boundary_x_packed_dev(hostptr , devptr, Iperiodic[0]);
	exchange_boundary_y_packed_dev(hostptr , devptr, Iperiodic[1]);
	exchange_boundary_z_packed_dev(hostptr , devptr, Iperiodic[2]);
}

void exchange_boundary_xyz_Async_packed_dev(REAL *hostptr , cudaField * devptr , hipStream_t *stream)
{
	exchange_boundary_x_Async_packed_dev(hostptr , devptr, Iperiodic[0], stream);
	exchange_boundary_y_Async_packed_dev(hostptr , devptr, Iperiodic[1], stream);
	exchange_boundary_z_Async_packed_dev(hostptr , devptr, Iperiodic[2], stream);
}


__global__ void cudaFieldBoundaryPack_kernel(cudaField data , cudaFieldPack pack , cudaJobPackage job){
    // eyes on cells WITH LAPs
	unsigned int x = blockDim.x * blockIdx.x + threadIdx.x;
	unsigned int y = blockDim.y * blockIdx.y + threadIdx.y;
    unsigned int z = blockDim.z * blockIdx.z + threadIdx.z;
	if(x < job.end.x && y < job.end.y && z < job.end.z){

        unsigned int pos = x + job.end.x*(y + job.end.y*z);
        x += job.start.x;
        y += job.start.y;
        z += job.start.z;
        *(pack.ptr + pos) = get_Field_LAP(data , x,y,z);
        
    }
}
__global__ void cudaFieldBoundaryUnpack_kernel(cudaField data , cudaFieldPack pack , cudaJobPackage job){
    // eyes on cells WITH LAPs
	unsigned int x = blockDim.x * blockIdx.x + threadIdx.x;
	unsigned int y = blockDim.y * blockIdx.y + threadIdx.y;
    unsigned int z = blockDim.z * blockIdx.z + threadIdx.z;
	if(x < job.end.x && y < job.end.y && z < job.end.z){

        unsigned int pos = x + job.end.x*(y + job.end.y*z);
        x += job.start.x;
        y += job.start.y;
        z += job.start.z;
        get_Field_LAP(data , x,y,z) = *(pack.ptr + pos);
        
    }
}


void cudaFieldBoundaryPack(cudaField * data , cudaFieldPack * pack, cudaJobPackage job_in)
{
    // job_in , to packed data , with LAP
    // job_in.start , job_in.size
    dim3 griddim , blockdim;
    cal_grid_block_dim(&griddim , &blockdim , BlockDimX , BlockDimY , BlockDimZ , job_in.end.x , job_in.end.y , job_in.end.z);
    CUDA_LAUNCH(( hipLaunchKernelGGL(cudaFieldBoundaryPack_kernel, dim3(griddim), dim3(blockdim), 0, 0, *data , *pack , job_in) ))
}

void cudaFieldBoundaryUnpack(cudaField * data , cudaFieldPack * pack , cudaJobPackage job_in){
    // job_in , to packed data , with LAP
    // job_in.start , job_in.size
    dim3 griddim , blockdim;
    cal_grid_block_dim(&griddim , &blockdim , BlockDimX , BlockDimY , BlockDimZ , job_in.end.x , job_in.end.y , job_in.end.z);
    CUDA_LAUNCH(( hipLaunchKernelGGL(cudaFieldBoundaryUnpack_kernel, dim3(griddim), dim3(blockdim), 0, 0, *data , *pack , job_in) ))
}

void exchange_boundary_x_packed_dev(REAL *hostptr , cudaField * devptr, int Iperiodic1)
{   
    cudaFieldPack * pack;
    MPI_Status status;
    int size = LAP*ny*nz;
    pack = b_xm;
    cudaJobPackage job(dim3(LAP,LAP,LAP),dim3(LAP,ny,nz));

    if(npx != 0 || Iperiodic1 == 1){
        cudaFieldBoundaryPack(devptr , pack ,job);
        CUDA_CALL(( hipMemcpy(pack_send_x , pack->ptr , size*sizeof(REAL) , hipMemcpyDeviceToHost) ))
    }
    MPI_Sendrecv(pack_send_x , size , OCFD_DATA_TYPE , ID_XM1 , 1 , pack_recv_x , size , OCFD_DATA_TYPE , ID_XP1 , 1 , MPI_COMM_WORLD , &status);
	if (npx != NPX0 - 1 || Iperiodic1 == 1){
        CUDA_CALL(( hipMemcpy(pack->ptr , pack_recv_x , size*sizeof(REAL) , hipMemcpyHostToDevice) ))
        job.start.x = nx_lap;
        cudaFieldBoundaryUnpack(devptr, pack ,job);
    }
    


    if(npx != NPX0 - 1 || Iperiodic1 == 1){
        job.start.x = nx;
        cudaFieldBoundaryPack(devptr , pack ,job);
        CUDA_CALL(( hipMemcpy(pack_send_x , pack->ptr , size*sizeof(REAL) , hipMemcpyDeviceToHost) ))
    }
    MPI_Sendrecv(pack_send_x , size , OCFD_DATA_TYPE , ID_XP1 , 1 , pack_recv_x , size , OCFD_DATA_TYPE , ID_XM1 , 1 , MPI_COMM_WORLD , &status);
	if (npx != 0 || Iperiodic1 == 1){
        CUDA_CALL(( hipMemcpy(pack->ptr , pack_recv_x , size*sizeof(REAL) , hipMemcpyHostToDevice) ))
        job.start.x = 0;
        cudaFieldBoundaryUnpack(devptr, pack ,job);
    }
    
}


void exchange_boundary_y_packed_dev(REAL *hostptr , cudaField * devptr, int Iperiodic1)
{   
    cudaFieldPack * pack;
    MPI_Status status;
    int size = LAP*nx*nz;
    pack = b_ym;
    cudaJobPackage job(dim3(LAP,LAP,LAP),dim3(nx , LAP ,nz));

    if(npy != 0 || Iperiodic1 == 1){
        cudaFieldBoundaryPack(devptr , pack ,job);
        CUDA_CALL(( hipMemcpy(pack_send_y , pack->ptr , size*sizeof(REAL) , hipMemcpyDeviceToHost) ))
    }
    MPI_Sendrecv(pack_send_y , size , OCFD_DATA_TYPE , ID_YM1 , 1 , pack_recv_y , size , OCFD_DATA_TYPE , ID_YP1 , 1 , MPI_COMM_WORLD , &status);
	if (npy != NPY0 - 1 || Iperiodic1 == 1){
        CUDA_CALL(( hipMemcpy(pack->ptr , pack_recv_y , size*sizeof(REAL) , hipMemcpyHostToDevice) ))
        job.start.y = ny_lap;
        cudaFieldBoundaryUnpack(devptr, pack ,job);
    }
    


    if(npy != NPY0 - 1 || Iperiodic1 == 1){
        job.start.y = ny;
        cudaFieldBoundaryPack(devptr , pack ,job);
        CUDA_CALL(( hipMemcpy(pack_send_y , pack->ptr , size*sizeof(REAL) , hipMemcpyDeviceToHost) ))
    }
    MPI_Sendrecv(pack_send_y , size , OCFD_DATA_TYPE , ID_YP1 , 1 , pack_recv_y , size , OCFD_DATA_TYPE , ID_YM1 , 1 , MPI_COMM_WORLD , &status);
	if (npy != 0 || Iperiodic1 == 1){
        CUDA_CALL(( hipMemcpy(pack->ptr , pack_recv_y , size*sizeof(REAL) , hipMemcpyHostToDevice) ))
        job.start.y = 0;
        cudaFieldBoundaryUnpack(devptr, pack ,job);
    }

}


void exchange_boundary_z_packed_dev(REAL *hostptr , cudaField * devptr, int Iperiodic1)
{   
    cudaFieldPack * pack;
    MPI_Status status;
    int size = LAP*nx*ny;
    pack = b_zm;
    cudaJobPackage job(dim3(LAP,LAP,LAP),dim3(nx,ny,LAP));

    if(npz != 0 || Iperiodic1 == 1){
        cudaFieldBoundaryPack(devptr , pack ,job);
        CUDA_CALL(( hipMemcpy(pack_send_z , pack->ptr , size*sizeof(REAL) , hipMemcpyDeviceToHost) ))
    }
    MPI_Sendrecv(pack_send_z , size , OCFD_DATA_TYPE , ID_ZM1 , 1 , pack_recv_z , size , OCFD_DATA_TYPE , ID_ZP1 , 1 , MPI_COMM_WORLD , &status);
	if (npz != NPZ0 - 1 || Iperiodic1 == 1){
        CUDA_CALL(( hipMemcpy(pack->ptr , pack_recv_z , size*sizeof(REAL) , hipMemcpyHostToDevice) ))
        job.start.z = nz_lap;
        cudaFieldBoundaryUnpack(devptr, pack ,job);
    }
    

    if(npz != NPZ0 - 1 || Iperiodic1 == 1){
        job.start.z = nz;
        cudaFieldBoundaryPack(devptr , pack ,job);
        CUDA_CALL(( hipMemcpy(pack_send_z , pack->ptr , size*sizeof(REAL) , hipMemcpyDeviceToHost) ))
    }
    MPI_Sendrecv(pack_send_z , size , OCFD_DATA_TYPE , ID_ZP1 , 1 , pack_recv_z , size , OCFD_DATA_TYPE , ID_ZM1 , 1 , MPI_COMM_WORLD , &status);
	if (npz != 0 || Iperiodic1 == 1){
        CUDA_CALL(( hipMemcpy(pack->ptr , pack_recv_z , size*sizeof(REAL) , hipMemcpyHostToDevice) ))
        job.start.z = 0;
        cudaFieldBoundaryUnpack(devptr, pack ,job);
    }

}

void cudaFieldBoundaryPack_Async(cudaField * data , cudaFieldPack * pack, cudaJobPackage job_in, hipStream_t *stream)
{
    // job_in , to packed data , with LAP
    // job_in.start , job_in.size
    dim3 griddim , blockdim;
    cal_grid_block_dim(&griddim , &blockdim , BlockDimX , BlockDimY , BlockDimZ , job_in.end.x , job_in.end.y , job_in.end.z);
    CUDA_LAUNCH(( hipLaunchKernelGGL(cudaFieldBoundaryPack_kernel, dim3(griddim), dim3(blockdim), 0, *stream, *data , *pack , job_in) ))
}

void cudaFieldBoundaryUnpack_Async(cudaField * data , cudaFieldPack * pack , cudaJobPackage job_in, hipStream_t *stream){
    // job_in , to packed data , with LAP
    // job_in.start , job_in.size
    dim3 griddim , blockdim;
    cal_grid_block_dim(&griddim , &blockdim , BlockDimX , BlockDimY , BlockDimZ , job_in.end.x , job_in.end.y , job_in.end.z);
    CUDA_LAUNCH(( hipLaunchKernelGGL(cudaFieldBoundaryUnpack_kernel, dim3(griddim), dim3(blockdim), 0, *stream, *data , *pack , job_in) ))
}

// 假设 , 仅仅交换边界
void exchange_boundary_x_Async_packed_dev(REAL *hostptr , cudaField * devptr, int Iperiodic1 , hipStream_t *stream)
{   
    cudaFieldPack * pack;
    MPI_Status status;
    int size = LAP*ny*nz;
    pack = b_xm;
    cudaJobPackage job(dim3(LAP,LAP,LAP),dim3(LAP,ny,nz));

    if(npx != 0 || Iperiodic1 == 1){
        cudaFieldBoundaryPack_Async(devptr , pack ,job, stream);
        CUDA_CALL(( hipMemcpyAsync(pack_send_x , pack->ptr , size*sizeof(REAL) , hipMemcpyDeviceToHost, *stream) ))
    }
    hipStreamSynchronize(*stream);
    MPI_Sendrecv(pack_send_x , size , OCFD_DATA_TYPE , ID_XM1 , 1 , pack_recv_x , size , OCFD_DATA_TYPE , ID_XP1 , 1 , MPI_COMM_WORLD , &status);
	if (npx != NPX0 - 1 || Iperiodic1 == 1){
        CUDA_CALL(( hipMemcpyAsync(pack->ptr , pack_recv_x , size*sizeof(REAL) , hipMemcpyHostToDevice, *stream) ))
        job.start.x = nx_lap;
        cudaFieldBoundaryUnpack_Async(devptr, pack ,job, stream);
    }
    


    if(npx != NPX0 - 1 || Iperiodic1 == 1){
        job.start.x = nx;
        cudaFieldBoundaryPack_Async(devptr , pack ,job, stream);
        CUDA_CALL(( hipMemcpyAsync(pack_send_x , pack->ptr , size*sizeof(REAL) , hipMemcpyDeviceToHost, *stream) ))
    }
    hipStreamSynchronize(*stream);
    MPI_Sendrecv(pack_send_x , size , OCFD_DATA_TYPE , ID_XP1 , 1 , pack_recv_x , size , OCFD_DATA_TYPE , ID_XM1 , 1 , MPI_COMM_WORLD , &status);
	if (npx != 0 || Iperiodic1 == 1){
        CUDA_CALL(( hipMemcpyAsync(pack->ptr , pack_recv_x , size*sizeof(REAL) , hipMemcpyHostToDevice, *stream) ))
        job.start.x = 0;
        cudaFieldBoundaryUnpack_Async(devptr, pack ,job, stream);
    }
    
}

void exchange_boundary_y_Async_packed_dev(REAL *hostptr , cudaField * devptr, int Iperiodic1 , hipStream_t *stream)
{   
    cudaFieldPack * pack;
    MPI_Status status;
    int size = LAP*nx*nz;
    pack = b_ym;
    cudaJobPackage job(dim3(LAP,LAP,LAP),dim3(nx , LAP ,nz));

    if(npy != 0 || Iperiodic1 == 1){
        cudaFieldBoundaryPack_Async(devptr , pack ,job, stream);
        CUDA_CALL(( hipMemcpyAsync(pack_send_y , pack->ptr , size*sizeof(REAL) , hipMemcpyDeviceToHost, *stream) ))
    }
    hipStreamSynchronize(*stream);
    MPI_Sendrecv(pack_send_y , size , OCFD_DATA_TYPE , ID_YM1 , 1 , pack_recv_y , size , OCFD_DATA_TYPE , ID_YP1 , 1 , MPI_COMM_WORLD , &status);
	if (npy != NPY0 - 1 || Iperiodic1 == 1){
        CUDA_CALL(( hipMemcpyAsync(pack->ptr , pack_recv_y , size*sizeof(REAL) , hipMemcpyHostToDevice, *stream) ))
        job.start.y = ny_lap;
        cudaFieldBoundaryUnpack_Async(devptr, pack ,job, stream);
    }
    


    if(npy != NPY0 - 1 || Iperiodic1 == 1){
        job.start.y = ny;
        cudaFieldBoundaryPack_Async(devptr , pack ,job, stream);
        CUDA_CALL(( hipMemcpyAsync(pack_send_y , pack->ptr , size*sizeof(REAL) , hipMemcpyDeviceToHost, *stream) ))
    }
    hipStreamSynchronize(*stream);
    MPI_Sendrecv(pack_send_y , size , OCFD_DATA_TYPE , ID_YP1 , 1 , pack_recv_y , size , OCFD_DATA_TYPE , ID_YM1 , 1 , MPI_COMM_WORLD , &status);
	if (npy != 0 || Iperiodic1 == 1){
        CUDA_CALL(( hipMemcpyAsync(pack->ptr , pack_recv_y , size*sizeof(REAL) , hipMemcpyHostToDevice, *stream) ))
        job.start.y = 0;
        cudaFieldBoundaryUnpack_Async(devptr, pack ,job, stream);
    }

}

void exchange_boundary_z_Async_packed_dev(REAL *hostptr , cudaField * devptr, int Iperiodic1 , hipStream_t *stream)
{   
    cudaFieldPack * pack;
    MPI_Status status;
    int size = LAP*nx*ny;
    pack = b_zm;
    cudaJobPackage job(dim3(LAP,LAP,LAP),dim3(nx,ny,LAP));

    if(npz != 0 || Iperiodic1 == 1){
        cudaFieldBoundaryPack_Async(devptr , pack ,job, stream);
        CUDA_CALL(( hipMemcpyAsync(pack_send_z , pack->ptr , size*sizeof(REAL) , hipMemcpyDeviceToHost, *stream) ))
    }
    hipStreamSynchronize(*stream);
    MPI_Sendrecv(pack_send_z , size , OCFD_DATA_TYPE , ID_ZM1 , 1 , pack_recv_z , size , OCFD_DATA_TYPE , ID_ZP1 , 1 , MPI_COMM_WORLD , &status);
	if (npz != NPZ0 - 1 || Iperiodic1 == 1){
        CUDA_CALL(( hipMemcpyAsync(pack->ptr , pack_recv_z , size*sizeof(REAL) , hipMemcpyHostToDevice, *stream) ))
        job.start.z = nz_lap;
        cudaFieldBoundaryUnpack_Async(devptr, pack ,job, stream);
    }
    


    if(npz != NPZ0 - 1 || Iperiodic1 == 1){
        job.start.z = nz;
        cudaFieldBoundaryPack_Async(devptr , pack ,job, stream);
        CUDA_CALL(( hipMemcpyAsync(pack_send_z , pack->ptr , size*sizeof(REAL) , hipMemcpyDeviceToHost, *stream) ))
    }
    hipStreamSynchronize(*stream);
    MPI_Sendrecv(pack_send_z , size , OCFD_DATA_TYPE , ID_ZP1 , 1 , pack_recv_z , size , OCFD_DATA_TYPE , ID_ZM1 , 1 , MPI_COMM_WORLD , &status);
	if (npz != 0 || Iperiodic1 == 1){
        CUDA_CALL(( hipMemcpyAsync(pack->ptr , pack_recv_z , size*sizeof(REAL) , hipMemcpyHostToDevice, *stream) ))
        job.start.z = 0;
        cudaFieldBoundaryUnpack_Async(devptr, pack ,job, stream);
    }

}


/* 

    switch(dir){
        case MPI_X_DIR : {
            MPI_Sendrecv(pack_send , size , OCFD_DATA_TYPE , ID_XM1 , 1 , pack_recv , size , OCFD_DATA_TYPE , ID_XP1 , 1 , MPI_COMM_WORKD , &status);
            break;
        }
        case MPI_Y_DIR : {
            MPI_Sendrecv(pack_send , size , OCFD_DATA_TYPE , ID_XM1 , 1 , pack_recv , size , OCFD_DATA_TYPE , ID_XP1 , 1 , MPI_COMM_WORKD , &status);
            break;
        }
        case MPI_Z_DIR : {
            MPI_Sendrecv(pack_send , size , OCFD_DATA_TYPE , ID_XM1 , 1 , pack_recv , size , OCFD_DATA_TYPE , ID_XP1 , 1 , MPI_COMM_WORKD , &status);
            break;
        }
    }

*/



#ifdef __cplusplus
}
#endif