1 /*
2 ! @file thesis_benchmark.c
3 ! @author Kyle Horne <horne.kyle@gmail.com>
4 ! @version 0.2
5 !
6 ! @section LICENSE
7 ! BSD style license
8 !
9 ! @section DESCRIPTION
10 ! Test program for pcgns library
11 */
12 
13 #include <stdio.h>
14 #include <string.h>
15 #include <stdlib.h>
16 #include <math.h>
17 
18 #include "pcgnslib.h"
19 #include "mpi.h"
20 
21 int comm_size;
22 int comm_rank;
23 MPI_Info info;
24 
25 double data_size = 1.0;
26 int N;
27 int Nl;
28 int pc;
29 int zpp = 2;
30 int ppz = 2;
31 int zc;
32 
33 int* zones;
34 int* subzones;
35 
36 double* x;
37 double* y;
38 double* z;
39 
40 cgsize_t* e;
41 
42 double* u;
43 double* v;
44 double* w;
45 
46 double* h;
47 
read_inputs(int * argc,char *** argv)48 int read_inputs(int* argc, char*** argv) {
49 	int k;
50 	if(comm_rank==0) {
51 		for(k=1;k<*argc;k++) {
52 			if(strcmp((*argv)[k],"-ds")==0) {
53 				k++;
54 				sscanf((*argv)[k],"%lf",&data_size);
55 				printf("data_size=%lf\n",data_size);
56 				}
57 			if(strcmp((*argv)[k],"-zpp")==0) {
58 				k++;
59 				sscanf((*argv)[k],"%d",&zpp);
60 				printf("zpp=%d\n",zpp);
61 				}
62 			if(strcmp((*argv)[k],"-ppz")==0) {
63 				k++;
64 				sscanf((*argv)[k],"%d",&ppz);
65 				printf("ppz=%d\n",ppz);
66 				}
67 			}
68 		}
69 	MPI_Bcast(&data_size,1,MPI_DOUBLE,0,MPI_COMM_WORLD);
70 	MPI_Bcast(&zpp,1,MPI_INT,0,MPI_COMM_WORLD);
71 	MPI_Bcast(&ppz,1,MPI_INT,0,MPI_COMM_WORLD);
72 	return 0;
73 	}
74 
initialize(int * argc,char *** argv)75 int initialize(int* argc, char*** argv) {
76 	int j,k;
77 	double theta;
78 	double r;
79 
80 	MPI_Init(argc,argv);
81 	MPI_Comm_size(MPI_COMM_WORLD, &comm_size);
82 	MPI_Comm_rank(MPI_COMM_WORLD, &comm_rank);
83 	MPI_Info_create(&info);
84 
85 	read_inputs(argc,argv);
86 
87 	N = (int)((data_size*1024*1024)/((double) sizeof(double)));
88 	pc = comm_size;
89 	zc = (pc*zpp)/ppz;
90 	Nl = N/comm_size/zpp;
91 
92 	zones = malloc(zc*sizeof(int));
93 	subzones = malloc(zpp*sizeof(int));
94 	for(k=0;k<zpp;k++) {
95 		zones[k] = comm_rank/ppz*zpp+k;
96 		subzones[k] = comm_rank%ppz;
97 		}
98 
99 	/* Initialize Arrays */
100 	x = malloc(Nl*sizeof(double));
101 	y = malloc(Nl*sizeof(double));
102 	z = malloc(Nl*sizeof(double));
103 
104 	e = malloc(Nl*sizeof(cgsize_t));
105 
106 	u = malloc(Nl*sizeof(double));
107 	v = malloc(Nl*sizeof(double));
108 	w = malloc(Nl*sizeof(double));
109 
110 	h = malloc(Nl*sizeof(double));
111 
112 	for(k=0;k<Nl;k++) {
113 		j = Nl*subzones[0]+k;
114 		theta = ((double) j)/((double) Nl*zpp);
115 		r = theta;
116 		x[k] = r*cos(theta);
117 		y[k] = r*sin(theta);
118 		z[k] = r;
119 		e[k] = j+1;
120 		u[k] = x[k];
121 		v[k] = y[k];
122 		w[k] = z[k];
123 		h[k] = r;
124 		}
125 
126 #if 0
127 	printf("%d: Nl %d\n",comm_rank,Nl);
128 	for(k=0;k<zpp;k++) printf("%d: Z%d.%d\n",comm_rank,zones[k],subzones[k]);
129 #endif
130 
131 	return 0;
132 	}
133 
finalize(void)134 int finalize(void) {
135 	free(zones);
136 	free(subzones);
137 
138 	free(x);
139 	free(y);
140 	free(z);
141 	free(e);
142 	free(u);
143 	free(v);
144 	free(w);
145 
146 	MPI_Finalize();
147 	return 0;
148 	}
149 
main(int argc,char * argv[])150 int main(int argc, char* argv[]) {
151 	int k;
152 	int F;
153 	int B;
154 
155 	int *Z, *E, *S, *A;
156 	int *Cx, *Cy, *Cz;
157 	int *Fu, *Fv, *Fw;
158 
159 	cgsize_t nijk[3];
160 	cgsize_t min, max;
161 
162 	double T0,T1;
163 	double Tw0,Tw1;
164 	double Tr0,Tr1;
165 	double t0,t1;
166 
167 	initialize(&argc,&argv);
168 
169 	Z = (int *)malloc(10*zc*sizeof(int));
170 	E = Z + zc;
171 	S = E + zc;
172 	A = S + zc;
173 	Cx = A + zc;
174 	Cy = Cx + zc;
175 	Cz = Cy + zc;
176 	Fu = Cz + zc;
177 	Fv = Fu + zc;
178 	Fw = Fv + zc;
179 
180 	if (cgp_open("thesis_benchmark.cgns", CG_MODE_WRITE, &F) ||
181 	    cg_base_write(F,"Base",3,3,&B))
182 	    cgp_error_exit();
183 
184 	nijk[0] = Nl*ppz;
185 	nijk[1] = Nl*ppz;
186 	nijk[2] = 0;
187 
188 	for(k=0;k<zc;k++) {
189 		char zonename[100+1];
190 		sprintf(zonename,"%s %d","Zone",k);
191 		if (cg_zone_write(F,B,zonename,nijk,CGNS_ENUMV(Unstructured),&(Z[k])) ||
192 		    cgp_coord_write(F,B,Z[k],CGNS_ENUMV(RealDouble),"CoordinateX",&(Cx[k])) ||
193 		    cgp_coord_write(F,B,Z[k],CGNS_ENUMV(RealDouble),"CoordinateY",&(Cy[k])) ||
194 		    cgp_coord_write(F,B,Z[k],CGNS_ENUMV(RealDouble),"CoordinateZ",&(Cz[k])) ||
195 		    cgp_section_write(F,B,Z[k],"Elements",CGNS_ENUMV(NODE),1,Nl*ppz,0,&(E[k])) ||
196 		    cg_sol_write(F,B,Z[k],"Solution",CGNS_ENUMV(Vertex),&S[k]) ||
197 		    cgp_field_write(F,B,Z[k],S[k],CGNS_ENUMV(RealDouble),"MomentumX",&(Fu[k])) ||
198 		    cgp_field_write(F,B,Z[k],S[k],CGNS_ENUMV(RealDouble),"MomentumY",&(Fv[k])) ||
199 		    cgp_field_write(F,B,Z[k],S[k],CGNS_ENUMV(RealDouble),"MomentumZ",&(Fw[k])))
200 		    cgp_error_exit();
201 		if (cg_goto(F,B,zonename,0,NULL) ||
202 		    cg_user_data_write("User Data") ||
203 		    cg_gorel(F, "User Data", 0, NULL) ||
204 		    cgp_array_write("phi",CGNS_ENUMV(RealDouble),1,nijk,&A[k]))
205 		    cgp_error_exit();
206 		}
207 
208 	MPI_Barrier(MPI_COMM_WORLD);
209 	T0 = MPI_Wtime();
210 	Tw0 = MPI_Wtime();
211 
212 	/* Writes */
213 	MPI_Barrier(MPI_COMM_WORLD);
214 	t0 = MPI_Wtime();
215 	for(k=0;k<zpp;k++) {
216 		min = subzones[k]*Nl+1;
217 		max = (subzones[k]+1)*Nl;
218 
219 		if (cgp_coord_write_data(F,B,Z[zones[k]],Cx[zones[k]],&min,&max,x) ||
220 		    cgp_coord_write_data(F,B,Z[zones[k]],Cy[zones[k]],&min,&max,y) ||
221 		    cgp_coord_write_data(F,B,Z[zones[k]],Cz[zones[k]],&min,&max,z))
222 		    cgp_error_exit();
223 		}
224 	MPI_Barrier(MPI_COMM_WORLD);
225 	t1 = MPI_Wtime();
226 	if(comm_rank==0) {
227 		printf("Coords Write\n");
228 		printf("\tTime=%lf\n",t1-t0);
229 		printf("\tBandwidth=%lf\n",3.0*data_size/(t1-t0));
230 		}
231 
232 	MPI_Barrier(MPI_COMM_WORLD);
233 	t0 = MPI_Wtime();
234 	for(k=0;k<zpp;k++) {
235 		min = subzones[k]*Nl+1;
236 		max = (subzones[k]+1)*Nl;
237 
238 		if (cgp_field_write_data(F,B,Z[zones[k]],S[zones[k]],Fu[zones[k]],&min,&max,u) ||
239 		    cgp_field_write_data(F,B,Z[zones[k]],S[zones[k]],Fv[zones[k]],&min,&max,v) ||
240 		    cgp_field_write_data(F,B,Z[zones[k]],S[zones[k]],Fw[zones[k]],&min,&max,w))
241 		    cgp_error_exit();
242 		}
243 	MPI_Barrier(MPI_COMM_WORLD);
244 	t1 = MPI_Wtime();
245 	if(comm_rank==0) {
246 		printf("Solutions Write\n");
247 		printf("\tTime=%lf\n",t1-t0);
248 		printf("\tBandwidth=%lf\n",3.0*data_size/(t1-t0));
249 		}
250 
251 	MPI_Barrier(MPI_COMM_WORLD);
252 	t0 = MPI_Wtime();
253 	for(k=0;k<zpp;k++) {
254 		min = subzones[k]*Nl+1;
255 		max = (subzones[k]+1)*Nl;
256 
257 		if (cg_goto(F,B,"Zone_t",Z[zones[k]],
258 		        "UserDefinedData_t",1,NULL) ||
259 		    cgp_array_write_data(A[zones[k]],&min,&max,h))
260 		    cgp_error_exit();
261 		}
262 	MPI_Barrier(MPI_COMM_WORLD);
263 	t1 = MPI_Wtime();
264 	if(comm_rank==0) {
265 		printf("Arrays Write\n");
266 		printf("\tTime=%lf\n",t1-t0);
267 		printf("\tBandwidth=%lf\n",data_size/(t1-t0));
268 		}
269 
270 	MPI_Barrier(MPI_COMM_WORLD);
271 	t0 = MPI_Wtime();
272 	for(k=0;k<zpp;k++) {
273 		min = subzones[k]*Nl+1;
274 		max = (subzones[k]+1)*Nl;
275 
276 		if (cgp_elements_write_data(F,B,Z[zones[k]],E[zones[k]],min,max,e))
277 		    cgp_error_exit();
278 		}
279 	MPI_Barrier(MPI_COMM_WORLD);
280 	t1 = MPI_Wtime();
281 	if(comm_rank==0) {
282 		printf("Elements Write\n");
283 		printf("\tTime=%lf\n",t1-t0);
284 		printf("\tBandwidth=%lf\n",((double) sizeof(int))/((double) sizeof(double))*data_size/(t1-t0));
285 		}
286 
287 	MPI_Barrier(MPI_COMM_WORLD);
288 	Tw1 = MPI_Wtime();
289 	if(comm_rank==0) {
290 		printf("Total Write Time=%lf\n",Tw1-Tw0);
291 		printf("Total Write Bandwidth=%lf\n",(6.0+((double) sizeof(int))/((double) sizeof(double)))*data_size/(Tw1-Tw0));
292 		}
293 
294 	MPI_Barrier(MPI_COMM_WORLD);
295 	t0 = MPI_Wtime();
296 	cgp_close(F);
297 	MPI_Barrier(MPI_COMM_WORLD);
298 	t1 = MPI_Wtime();
299 	if(comm_rank==0) printf("Close_Time=%lf\n",t1-t0);
300 
301 	/*=======
302 	 *=Reads=
303 	 *=======*/
304 
305 	if (cgp_open("thesis_benchmark.cgns", CG_MODE_READ, &F))
306 	    cgp_error_exit();
307 
308 	MPI_Barrier(MPI_COMM_WORLD);
309 	Tr0 = MPI_Wtime();
310 
311 	MPI_Barrier(MPI_COMM_WORLD);
312 	t0 = MPI_Wtime();
313 	for(k=0;k<zpp;k++) {
314 		min = subzones[k]*Nl+1;
315 		max = (subzones[k]+1)*Nl;
316 
317 		if (cgp_coord_read_data(F,B,Z[zones[k]],Cx[zones[k]],&min,&max,x) ||
318 		    cgp_coord_read_data(F,B,Z[zones[k]],Cy[zones[k]],&min,&max,y) ||
319 		    cgp_coord_read_data(F,B,Z[zones[k]],Cz[zones[k]],&min,&max,z))
320 		    cgp_error_exit();
321 		}
322 	MPI_Barrier(MPI_COMM_WORLD);
323 	t1 = MPI_Wtime();
324 	if(comm_rank==0) {
325 		printf("Coords Read\n");
326 		printf("\tTime=%lf\n",t1-t0);
327 		printf("\tBandwidth=%lf\n",3.0*data_size/(t1-t0));
328 		}
329 
330 	MPI_Barrier(MPI_COMM_WORLD);
331 	t0 = MPI_Wtime();
332 	for(k=0;k<zpp;k++) {
333 		min = subzones[k]*Nl+1;
334 		max = (subzones[k]+1)*Nl;
335 
336 		if (cgp_field_read_data(F,B,Z[zones[k]],S[zones[k]],Fu[zones[k]],&min,&max,u) ||
337 		    cgp_field_read_data(F,B,Z[zones[k]],S[zones[k]],Fv[zones[k]],&min,&max,v) ||
338 		    cgp_field_read_data(F,B,Z[zones[k]],S[zones[k]],Fw[zones[k]],&min,&max,w))
339 		    cgp_error_exit();
340 		}
341 	MPI_Barrier(MPI_COMM_WORLD);
342 	t1 = MPI_Wtime();
343 	if(comm_rank==0) {
344 		printf("Solutions Read\n");
345 		printf("\tTime=%lf\n",t1-t0);
346 		printf("\tBandwidth=%lf\n",3.0*data_size/(t1-t0));
347 		}
348 
349 	MPI_Barrier(MPI_COMM_WORLD);
350 	t0 = MPI_Wtime();
351 	for(k=0;k<zpp;k++) {
352 		min = subzones[k]*Nl+1;
353 		max = (subzones[k]+1)*Nl;
354 
355 		if (cg_goto(F,B,"Zone_t",Z[zones[k]],
356 		        "UserDefinedData_t",1,NULL) ||
357 		    cgp_array_read_data(A[zones[k]],&min,&max,h))
358 		    cgp_error_exit();
359 		}
360 	MPI_Barrier(MPI_COMM_WORLD);
361 	t1 = MPI_Wtime();
362 	if(comm_rank==0) {
363 		printf("Arrays Read\n");
364 		printf("\tTime=%lf\n",t1-t0);
365 		printf("\tBandwidth=%lf\n",data_size/(t1-t0));
366 		}
367 
368 	MPI_Barrier(MPI_COMM_WORLD);
369 	t0 = MPI_Wtime();
370 	for(k=0;k<zpp;k++) {
371 		min = subzones[k]*Nl+1;
372 		max = (subzones[k]+1)*Nl;
373 
374 		if (cgp_elements_read_data(F,B,Z[zones[k]],E[zones[k]],min,max,e))
375 		    cgp_error_exit();
376 		}
377 	MPI_Barrier(MPI_COMM_WORLD);
378 	t1 = MPI_Wtime();
379 	if(comm_rank==0) {
380 		printf("Elements Read\n");
381 		printf("\tTime=%lf\n",t1-t0);
382 		printf("\tBandwidth=%lf\n",((double) sizeof(int))/((double) sizeof(double))*data_size/(t1-t0));
383 		}
384 
385 	MPI_Barrier(MPI_COMM_WORLD);
386 	Tr1 = MPI_Wtime();
387 	if(comm_rank==0) {
388 		printf("Total Read Time=%lf\n",Tr1-Tr0);
389 		printf("Total Read Bandwidth=%lf\n",(6.0+((double) sizeof(int))/((double) sizeof(double)))*data_size/(Tr1-Tr0));
390 		}
391 
392 	MPI_Barrier(MPI_COMM_WORLD);
393 	t0 = MPI_Wtime();
394 	cgp_close(F);
395 	MPI_Barrier(MPI_COMM_WORLD);
396 	t1 = MPI_Wtime();
397 	if(comm_rank==0) printf("Close_Time=%lf\n",t1-t0);
398 
399 	MPI_Barrier(MPI_COMM_WORLD);
400 	T1 = MPI_Wtime();
401 
402 	if(comm_rank==0) {
403 		printf("Total Time=%lf\n",T1-T0);
404 		}
405 
406 	finalize();
407 
408 	return 0;
409 	}
410