-
Notifications
You must be signed in to change notification settings - Fork 31
/
Copy pathTestForceMoments.cpp
331 lines (281 loc) · 9.66 KB
/
TestForceMoments.cpp
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
//*************************************************************************
// Lattice Boltzmann Simulator for Single Phase Flow in Porous Media
// James E. McCLure
//*************************************************************************
#include <stdio.h>
#include <iostream>
#include <fstream>
#include "common/ScaLBL.h"
#include "common/MPI.h"
using namespace std;
extern void PrintNeighborList(int * neighborList, int Np, int rank) {
if (rank == 0) {
int neighbor;
for (int i = 0; i < Np; i++) {
printf("idx=%d: ",i);
for (int l = 0; l < 10; l++) { // was 18
neighbor = neighborList[l*Np + i];
printf("%d ",neighbor);
}
printf("\n");
}
printf("\n\n");
}
}
std::shared_ptr<Database> loadInputs( int nprocs )
{
auto db = std::make_shared<Database>( "Domain.in" );
db->putScalar<int>( "BC", 0 );
if ( nprocs == 1 ){
db->putVector<int>( "nproc", { 1, 1, 1 } );
db->putVector<int>( "n", { 3, 3, 3 } );
db->putScalar<int>( "nspheres", 0 );
db->putVector<double>( "L", { 1, 1, 1 } );
}
return db;
}
//***************************************************************************************
int main(int argc, char **argv)
{
// Initialize MPI
Utilities::startup( argc, argv );
Utilities::MPI comm( MPI_COMM_WORLD );
int rank = comm.getRank();
int nprocs = comm.getSize();
int check=0;
{
// parallel domain size (# of sub-domains)
int iproc,jproc,kproc;
if (rank == 0){
printf("********************************************************\n");
printf("Running Unit Test: TestForceMoments \n");
printf("********************************************************\n");
}
// BGK Model parameters
double tau,Fx,Fy,Fz;
// Domain variables
int i,j,k,n;
int dim = 3; if (rank == 0) printf("dim=%d\n",dim);
int timestep = 0;
int timesteps = 2;
tau =1.0;
double rlx_setA=1.0/tau;
double rlx_setB = 8.f*(2.f-rlx_setA)/(8.f-rlx_setA);
Fx = Fy = 1.0;
Fz = 1.0;
// Load inputs
string FILENAME = argv[1];
// Load inputs
if (rank==0) printf("Loading input database \n");
auto db = std::make_shared<Database>( FILENAME );
auto domain_db = db->getDatabase( "Domain" );
int Nx = domain_db->getVector<int>( "n" )[0];
int Ny = domain_db->getVector<int>( "n" )[1];
int Nz = domain_db->getVector<int>( "n" )[2];
int nprocx = domain_db->getVector<int>( "nproc" )[0];
int nprocy = domain_db->getVector<int>( "nproc" )[1];
int nprocz = domain_db->getVector<int>( "nproc" )[2];
if (rank==0){
printf("********************************************************\n");
printf("Sub-domain size = %i x %i x %i\n",Nx,Ny,Nz);
printf("********************************************************\n");
}
comm.barrier();
kproc = rank/(nprocx*nprocy);
jproc = (rank-nprocx*nprocy*kproc)/nprocx;
iproc = rank-nprocx*nprocy*kproc-nprocz*jproc;
if (rank == 0) {
printf("i,j,k proc=%d %d %d \n",iproc,jproc,kproc);
}
comm.barrier();
if (rank == 1){
printf("i,j,k proc=%d %d %d \n",iproc,jproc,kproc);
printf("\n\n");
}
std::shared_ptr<Domain> Dm(new Domain(domain_db,comm));
Nx += 2;
Ny += 2;
Nz += 2;
int N = Nx*Ny*Nz;
//.......................................................................
// Assign the phase ID field
//.......................................................................
char LocalRankString[8];
sprintf(LocalRankString,"%05d",rank);
char LocalRankFilename[40];
sprintf(LocalRankFilename,"ID.%05i",rank);
/*
FILE *IDFILE = fopen(LocalRankFilename,"rb");
if (IDFILE==NULL) ERROR("Error opening file: ID.xxxxx");
fread(Dm->id,1,N,IDFILE);
fclose(IDFILE);
*/
// initialize empty domain
for (k=0;k<Nz;k++){
for (j=0;j<Ny;j++){
for (i=0;i<Nx;i++){
n = k*Nx*Ny+j*Nx+i;
Dm->id[n]=1;
}
}
}
Dm->CommInit();
comm.barrier();
if (rank == 0) cout << "Domain set." << endl;
int Np=0; // number of local pore nodes
for (k=1;k<Nz-1;k++){
for (j=1;j<Ny-1;j++){
for (i=1;i<Nx-1;i++){
n = k*Nx*Ny+j*Nx+i;
if (Dm->id[n] > 0){
Np++;
}
}
}
}
if (rank==0) printf ("Create ScaLBL_Communicator \n");
auto ScaLBL_Comm = std::shared_ptr<ScaLBL_Communicator>(new ScaLBL_Communicator(Dm));
//...........device phase ID.................................................
if (rank==0) printf ("Copying phase ID to device \n");
char *ID;
ScaLBL_AllocateDeviceMemory((void **) &ID, N); // Allocate device memory
// Copy to the device
ScaLBL_CopyToDevice(ID, Dm->id.data(), N);
//...........................................................................
if (rank==0){
printf("Total domain size = %i \n",N);
printf("Reduced domain size = %i \n",Np);
}
// LBM variables
if (rank==0) printf ("Allocating distributions \n");
int neighborSize=18*Np*sizeof(int);
int *neighborList;
IntArray Map(Nx,Ny,Nz);
neighborList= new int[18*Np];
ScaLBL_Comm->MemoryOptimizedLayoutAA(Map,neighborList,Dm->id.data(),Np,1);
if (rank == 0) PrintNeighborList(neighborList,Np, rank);
comm.barrier();
//......................device distributions.................................
int dist_mem_size = Np*sizeof(double);
int *NeighborList;
// double *f_even,*f_odd;
double * dist;
double * Velocity;
//...........................................................................
ScaLBL_AllocateDeviceMemory((void **) &dist, 19*dist_mem_size);
ScaLBL_AllocateDeviceMemory((void **) &NeighborList, neighborSize);
ScaLBL_AllocateDeviceMemory((void **) &Velocity, 3*sizeof(double)*Np);
ScaLBL_CopyToDevice(NeighborList, neighborList, neighborSize);
//...........................................................................
/*
* AA Algorithm begins here
*
*/
ScaLBL_D3Q19_Init(dist, Np);
//.......create and start timer............
double starttime,stoptime,cputime;
ScaLBL_DeviceBarrier(); comm.barrier();
starttime = Utilities::MPI::time();
/************ MAIN ITERATION LOOP (timing communications)***************************************/
//ScaLBL_Comm->SendD3Q19(dist, &dist[10*Np]);
//ScaLBL_Comm->RecvD3Q19(dist, &dist[10*Np]);
ScaLBL_DeviceBarrier(); comm.barrier();
if (rank==0) printf("Beginning AA timesteps...\n");
if (rank==0) printf("********************************************************\n");
if (rank==0) printf("No. of timesteps for timing: %i \n", timesteps);
while (timestep < 2) {
ScaLBL_Comm->SendD3Q19AA(dist); //READ FROM NORMAL
ScaLBL_D3Q19_AAodd_MRT(NeighborList, dist, ScaLBL_Comm->first_interior, ScaLBL_Comm->last_interior, Np, rlx_setA, rlx_setB, Fx, Fy, Fz);
ScaLBL_Comm->RecvD3Q19AA(dist); //WRITE INTO OPPOSITE
ScaLBL_D3Q19_AAodd_MRT(NeighborList, dist, 0, ScaLBL_Comm->next, Np, rlx_setA, rlx_setB, Fx, Fy, Fz);
ScaLBL_DeviceBarrier(); comm.barrier();
timestep++;
ScaLBL_Comm->SendD3Q19AA(dist); //READ FORM NORMAL
ScaLBL_D3Q19_AAeven_MRT(dist, ScaLBL_Comm->first_interior, ScaLBL_Comm->last_interior, Np, rlx_setA, rlx_setB, Fx, Fy, Fz);
ScaLBL_Comm->RecvD3Q19AA(dist); //WRITE INTO OPPOSITE
ScaLBL_D3Q19_AAeven_MRT(dist, 0, ScaLBL_Comm->next, Np, rlx_setA, rlx_setB, Fx, Fy, Fz);
ScaLBL_DeviceBarrier(); comm.barrier();
timestep++;
//************************************************************************/
timestep++;
}
//************************************************************************/
stoptime = Utilities::MPI::time();
// cout << "CPU time: " << (stoptime - starttime) << " seconds" << endl;
cputime = stoptime - starttime;
// cout << "Lattice update rate: "<< double(Nx*Ny*Nz*timestep)/cputime/1000000 << " MLUPS" << endl;
double MLUPS = double(Np*timestep)/cputime/1000000;
if (rank==0) printf("********************************************************\n");
if (rank==0) printf("CPU time = %f \n", cputime);
if (rank==0) printf("Lattice update rate (per process)= %f MLUPS \n", MLUPS);
MLUPS *= nprocs;
if (rank==0) printf("Lattice update rate (process)= %f MLUPS \n", MLUPS);
if (rank==0) printf("********************************************************\n");
int SIZE=Np*sizeof(double);
double *DIST;
DIST= new double [19*Np];
ScaLBL_CopyToHost(&DIST[0],&dist[0],19*SIZE);
i=Nx/2;
printf("x = constant \n");
for (int q=0; q<9; q++){
int a = 2*q+1;
int b = 2*(q+1);
printf("************* \n");
printf("print slice for distribution pair %i,%i \n",a,b);
for (k=1;k<Nz-1;k++){
for (j=1;j<Ny-1;j++){
n = k*Nx*Ny+j*Nx+i;
//printf("%i ",Dm->id[n]);
n = Map(i,j,k);
double fa = DIST[(2*q+1)*Np+n];
double fb = DIST[2*(q+1)*Np+n];
printf("%f,%f ",fa,fb);
}
printf("\n");
}
printf("************* \n");
}
printf("y = constant \n");
j=Ny/2;
for (int q=0; q<9; q++){
int a = 2*q+1;
int b = 2*(q+1);
printf("************* \n");
printf("print slice for distribution pair %i,%i \n",a,b);
for (k=1;k<Nz-1;k++){
for (i=1;i<Nx-1;i++){
n = k*Nx*Ny+j*Nx+i;
//printf("%i ",Dm->id[n]);
n = Map(i,j,k);
double fa = DIST[(2*q+1)*Np+n];
double fb = DIST[2*(q+1)*Np+n];
printf("%f,%f ",fa,fb);
}
printf("\n");
}
printf("************* \n");
}
k=Nz/2;
printf("z = constant \n");
for (int q=0; q<9; q++){
int a = 2*q+1;
int b = 2*(q+1);
printf("************* \n");
printf("print slice for distribution pair %i,%i \n",a,b);
for (j=1;j<Ny-1;j++){
for (i=1;i<Nx-1;i++){
n = k*Nx*Ny+j*Nx+i;
//printf("%i ",Dm->id[n]);
n = Map(i,j,k);
double fa = DIST[(2*q+1)*Np+n];
double fb = DIST[2*(q+1)*Np+n];
printf("%f,%f ",fa,fb);
}
printf("\n");
}
printf("************* \n");
}
}
Utilities::shutdown();
return check;
}