xref: /phasta/checkphasta/checkphasta.cpp (revision 16223cb9c3f88b34f2cb94151b5cf5ffc1aac5e2)
1*16223cb9SCameron Smith #include <cstdio>
2*16223cb9SCameron Smith #include <cstdlib>
3*16223cb9SCameron Smith #include <cassert>
4*16223cb9SCameron Smith #include <cstring>
5*16223cb9SCameron Smith #include <cmath>
6*16223cb9SCameron Smith #include <string>
7*16223cb9SCameron Smith #include <vector>
8*16223cb9SCameron Smith #include <algorithm>
9*16223cb9SCameron Smith #include <set>
10*16223cb9SCameron Smith 
11*16223cb9SCameron Smith #include <sys/types.h>
12*16223cb9SCameron Smith #include <sys/stat.h>
13*16223cb9SCameron Smith #include <dirent.h>
14*16223cb9SCameron Smith 
15*16223cb9SCameron Smith #define OMPI_SKIP_MPICXX 1
16*16223cb9SCameron Smith #include <mpi.h>
17*16223cb9SCameron Smith 
18*16223cb9SCameron Smith #include "phastaIO.h"
19*16223cb9SCameron Smith 
20*16223cb9SCameron Smith //Provided by phastaIO
21*16223cb9SCameron Smith void Gather_Headers( int* fileDescriptor, std::vector< std::string >& headers );
22*16223cb9SCameron Smith 
23*16223cb9SCameron Smith char read_solution(double** solutiono, int* size, int* nshgo, int* ndofo,
24*16223cb9SCameron Smith 		int nump, int rank, int timestep, char* casedir);
25*16223cb9SCameron Smith 
26*16223cb9SCameron Smith std::set<int>* find_timesteps(char* casedir, int nump);
27*16223cb9SCameron Smith double compare_solution(char* lpath, char* rpath, int timestep, int nump);
28*16223cb9SCameron Smith 
29*16223cb9SCameron Smith int main(int argc, char** argv)
30*16223cb9SCameron Smith {
31*16223cb9SCameron Smith 	int rank;
32*16223cb9SCameron Smith 	int size;
33*16223cb9SCameron Smith 	MPI_Init(&argc, &argv);
34*16223cb9SCameron Smith 
35*16223cb9SCameron Smith 	MPI_Comm_rank(MPI_COMM_WORLD, &rank);
36*16223cb9SCameron Smith 	MPI_Comm_size(MPI_COMM_WORLD, &size);
37*16223cb9SCameron Smith 
38*16223cb9SCameron Smith 	assert(argc>2);
39*16223cb9SCameron Smith 	char* lpath = argv[1];
40*16223cb9SCameron Smith 	char* rpath = argv[2];
41*16223cb9SCameron Smith 
42*16223cb9SCameron Smith 	int ndof;
43*16223cb9SCameron Smith 	int nshg;
44*16223cb9SCameron Smith 	int solsize;
45*16223cb9SCameron Smith 	double* solution;
46*16223cb9SCameron Smith 
47*16223cb9SCameron Smith 	std::set<int>* l_timesteps = find_timesteps(lpath, size);
48*16223cb9SCameron Smith 	std::set<int>* r_timesteps = find_timesteps(rpath, size);
49*16223cb9SCameron Smith 	std::set<int>* timesteps_to_check = new std::set<int>;
50*16223cb9SCameron Smith 	std::set_intersection(l_timesteps->begin(), l_timesteps->end(),
51*16223cb9SCameron Smith 			r_timesteps->begin(), r_timesteps->end(),
52*16223cb9SCameron Smith 			std::inserter(*timesteps_to_check, timesteps_to_check->begin()));
53*16223cb9SCameron Smith         delete l_timesteps;
54*16223cb9SCameron Smith         delete r_timesteps;
55*16223cb9SCameron Smith 	if(rank == 0)
56*16223cb9SCameron Smith 		printf("Found %d common timesteps\n",
57*16223cb9SCameron Smith 			       	timesteps_to_check->size());
58*16223cb9SCameron Smith #ifdef DBGONLY
59*16223cb9SCameron Smith 	read_solution(&solution, &solsize, &nshg, &ndof, size, rank, 0, "./");
60*16223cb9SCameron Smith 	printf("nshg: %d, ndof: %d\n", nshg, ndof);
61*16223cb9SCameron Smith 	assert(solsize == ndof*nshg);
62*16223cb9SCameron Smith #endif
63*16223cb9SCameron Smith 	double maxerror = 0.0;
64*16223cb9SCameron Smith 	double error;
65*16223cb9SCameron Smith 	double gblmaxerror;
66*16223cb9SCameron Smith 	for(std::set<int>::iterator i = timesteps_to_check->begin();
67*16223cb9SCameron Smith 			i!=timesteps_to_check->end();i++)
68*16223cb9SCameron Smith 	{
69*16223cb9SCameron Smith 		error = compare_solution(lpath, rpath, *i, size);
70*16223cb9SCameron Smith 		if(error>maxerror) maxerror = error;
71*16223cb9SCameron Smith 	}
72*16223cb9SCameron Smith         delete timesteps_to_check;
73*16223cb9SCameron Smith 	MPI_Barrier(MPI_COMM_WORLD);
74*16223cb9SCameron Smith 	MPI_Reduce(&maxerror, &gblmaxerror, 1, MPI_DOUBLE, MPI_MAX, 0,
75*16223cb9SCameron Smith 		MPI_COMM_WORLD);
76*16223cb9SCameron Smith 	if(rank == 0) printf("Maximum difference across all timesteps: %e\n",
77*16223cb9SCameron Smith 			gblmaxerror);
78*16223cb9SCameron Smith 	MPI_Finalize();
79*16223cb9SCameron Smith }
80*16223cb9SCameron Smith double compare_solution(char* lpath, char* rpath, int timestep, int nump)
81*16223cb9SCameron Smith {
82*16223cb9SCameron Smith 	int rank;
83*16223cb9SCameron Smith 	MPI_Comm_rank(MPI_COMM_WORLD, &rank);
84*16223cb9SCameron Smith 	double* lsol;
85*16223cb9SCameron Smith 	double* rsol;
86*16223cb9SCameron Smith 	int lsize;
87*16223cb9SCameron Smith 	int rsize;
88*16223cb9SCameron Smith 
89*16223cb9SCameron Smith 	read_solution(&lsol, &lsize, NULL, NULL, nump, rank, timestep, lpath);
90*16223cb9SCameron Smith 	read_solution(&rsol, &rsize, NULL, NULL, nump, rank, timestep, rpath);
91*16223cb9SCameron Smith 
92*16223cb9SCameron Smith 	double maxdiff=0.0;
93*16223cb9SCameron Smith 	double gblmaxdiff;
94*16223cb9SCameron Smith 	if(lsize != rsize)
95*16223cb9SCameron Smith 	{
96*16223cb9SCameron Smith 		printf("Error: Solution sizes different: %d, %d\n",
97*16223cb9SCameron Smith 				lsize, rsize);
98*16223cb9SCameron Smith 		assert(lsize == rsize);
99*16223cb9SCameron Smith 	}
100*16223cb9SCameron Smith 	for(int i=0;i<lsize;i++)
101*16223cb9SCameron Smith 	{
102*16223cb9SCameron Smith 		double diff = fabs(lsol[i]-rsol[i]);
103*16223cb9SCameron Smith 		if(diff > maxdiff) maxdiff = diff;
104*16223cb9SCameron Smith 	}
105*16223cb9SCameron Smith         free(lsol);
106*16223cb9SCameron Smith         free(rsol);
107*16223cb9SCameron Smith 	MPI_Reduce(&maxdiff, &gblmaxdiff, 1, MPI_DOUBLE, MPI_MAX, 0,
108*16223cb9SCameron Smith 		       	MPI_COMM_WORLD);
109*16223cb9SCameron Smith 	MPI_Barrier(MPI_COMM_WORLD); //TODO: debugging only
110*16223cb9SCameron Smith 	if(rank == 0)
111*16223cb9SCameron Smith 	printf("Timestep: %d, maximum difference: %e\n", timestep, gblmaxdiff);
112*16223cb9SCameron Smith 	return gblmaxdiff;
113*16223cb9SCameron Smith 
114*16223cb9SCameron Smith }
115*16223cb9SCameron Smith char read_solution(double** solutiono, int* size, int* nshgo, int* ndofo,
116*16223cb9SCameron Smith 		int nump, int rank, int timestep, char* casedir)
117*16223cb9SCameron Smith {
118*16223cb9SCameron Smith 	int iarray[10];
119*16223cb9SCameron Smith         const char* iformat = "binary";
120*16223cb9SCameron Smith         int ithree=3;
121*16223cb9SCameron Smith         int igeombc;
122*16223cb9SCameron Smith         char* fn;
123*16223cb9SCameron Smith 	int nshg;
124*16223cb9SCameron Smith 	int ndof;
125*16223cb9SCameron Smith 	double* solution;
126*16223cb9SCameron Smith 	if(nump == 1)
127*16223cb9SCameron Smith 		asprintf(&fn, "%s/restart.%d.%d",
128*16223cb9SCameron Smith 				casedir,timestep,rank+1);
129*16223cb9SCameron Smith 	else
130*16223cb9SCameron Smith 		asprintf(&fn, "%s/%d-procs_case/restart.%d.%d",
131*16223cb9SCameron Smith 				casedir, nump,timestep,rank+1);
132*16223cb9SCameron Smith         openfile(fn, "read", &igeombc);
133*16223cb9SCameron Smith 	//TODO: error handle
134*16223cb9SCameron Smith 	readheader(&igeombc, "solution", (void*) iarray, &ithree, "integer", iformat);
135*16223cb9SCameron Smith 	nshg = iarray[0];
136*16223cb9SCameron Smith 	ndof = iarray[1];
137*16223cb9SCameron Smith 	if(size != NULL)
138*16223cb9SCameron Smith 		*size = nshg*ndof;
139*16223cb9SCameron Smith 	solution = (double*) malloc(sizeof(double)*nshg*ndof);
140*16223cb9SCameron Smith 	readdatablock(&igeombc, "solution", solution, size, "double", iformat);
141*16223cb9SCameron Smith 	closefile(&igeombc, "read");
142*16223cb9SCameron Smith 	if(solutiono != NULL)
143*16223cb9SCameron Smith 		*solutiono = solution;
144*16223cb9SCameron Smith 	if(nshgo != NULL)
145*16223cb9SCameron Smith 		*nshgo = nshg;
146*16223cb9SCameron Smith 	if(ndofo != NULL)
147*16223cb9SCameron Smith 		*ndofo = ndof;
148*16223cb9SCameron Smith 	free(fn);
149*16223cb9SCameron Smith 	return(0);
150*16223cb9SCameron Smith }
151*16223cb9SCameron Smith 
152*16223cb9SCameron Smith std::set<int>* find_timesteps(char* casedir, int nump)
153*16223cb9SCameron Smith {
154*16223cb9SCameron Smith 	char* path;
155*16223cb9SCameron Smith 	char* fullpath;
156*16223cb9SCameron Smith 	DIR* dir;
157*16223cb9SCameron Smith 	struct dirent* d;
158*16223cb9SCameron Smith 	int part, ts;
159*16223cb9SCameron Smith 	std::set<int>* step_list = new std::set<int>;
160*16223cb9SCameron Smith 
161*16223cb9SCameron Smith 	if(nump == 1)
162*16223cb9SCameron Smith 		asprintf(&path, "%s", casedir);
163*16223cb9SCameron Smith 	else
164*16223cb9SCameron Smith 		asprintf(&path, "%s/%d-procs_case", casedir, nump);
165*16223cb9SCameron Smith 	dir = opendir(path);
166*16223cb9SCameron Smith 	if(!dir)
167*16223cb9SCameron Smith 	{
168*16223cb9SCameron Smith 		perror("Error opening case: ");
169*16223cb9SCameron Smith 		MPI_Abort(MPI_COMM_WORLD,1);
170*16223cb9SCameron Smith 	}
171*16223cb9SCameron Smith 	while((d=readdir(dir)))
172*16223cb9SCameron Smith 	{
173*16223cb9SCameron Smith 		asprintf(&fullpath, "%s/%s", path, d->d_name);
174*16223cb9SCameron Smith 		if(sscanf(d->d_name, "restart.%d.%d", &ts, &part)==2)
175*16223cb9SCameron Smith 		{
176*16223cb9SCameron Smith 			step_list->insert(ts);
177*16223cb9SCameron Smith 		}
178*16223cb9SCameron Smith 	}
179*16223cb9SCameron Smith 	return(step_list);
180*16223cb9SCameron Smith 	free(path);
181*16223cb9SCameron Smith }
182*16223cb9SCameron Smith 
183*16223cb9SCameron Smith void check_ilwork(int* ilwork, int nlwork)
184*16223cb9SCameron Smith {
185*16223cb9SCameron Smith 	int numtask = ilwork[0];
186*16223cb9SCameron Smith 	int itkbeg = 0; //task offset
187*16223cb9SCameron Smith 	printf("%d tasks\n", numtask);
188*16223cb9SCameron Smith 
189*16223cb9SCameron Smith 	for(int i=0;i<numtask;i++)
190*16223cb9SCameron Smith 	{
191*16223cb9SCameron Smith 		int itag = ilwork[itkbeg+1]; //mpi tag
192*16223cb9SCameron Smith 		int iacc = ilwork[itkbeg+2]; //0 for slave, 1 for master
193*16223cb9SCameron Smith 		int iother = ilwork[itkbeg+3]-1; //other rank (see ctypes.f for off by one)
194*16223cb9SCameron Smith 		int numseg = ilwork[itkbeg+4]; //number of segments
195*16223cb9SCameron Smith 		printf("Comm with rank: %d\n", iother);
196*16223cb9SCameron Smith 		for(int j=0;j<numseg;j++)
197*16223cb9SCameron Smith 		{
198*16223cb9SCameron Smith 			int isgbeg = ilwork[itkbeg+5+(j*2)]; //first idx of seg
199*16223cb9SCameron Smith 			int lenseg = ilwork[itkbeg+6+(j*2)]; //length of seg
200*16223cb9SCameron Smith 
201*16223cb9SCameron Smith 			printf("isgbeg: %d, len: %d\n", isgbeg, lenseg);
202*16223cb9SCameron Smith 			assert(itkbeg+6+(j*2) < nlwork);
203*16223cb9SCameron Smith 		}
204*16223cb9SCameron Smith 		itkbeg+= 4+2*numseg;
205*16223cb9SCameron Smith 	}
206*16223cb9SCameron Smith }
207*16223cb9SCameron Smith 
208*16223cb9SCameron Smith void read_ilwork()
209*16223cb9SCameron Smith {
210*16223cb9SCameron Smith 	int rank;
211*16223cb9SCameron Smith 	int size;
212*16223cb9SCameron Smith 
213*16223cb9SCameron Smith 	MPI_Comm_rank(MPI_COMM_WORLD, &rank);
214*16223cb9SCameron Smith 	MPI_Comm_size(MPI_COMM_WORLD, &size);
215*16223cb9SCameron Smith 
216*16223cb9SCameron Smith 	int iarray[10];
217*16223cb9SCameron Smith 	const char* iformat = "binary";
218*16223cb9SCameron Smith 	int ione=1;
219*16223cb9SCameron Smith 	int itwo=2;
220*16223cb9SCameron Smith 	int igeombc;
221*16223cb9SCameron Smith 	int nilwork=0;
222*16223cb9SCameron Smith 	int tmp;
223*16223cb9SCameron Smith 	int* ilwork;
224*16223cb9SCameron Smith 	char* fn;
225*16223cb9SCameron Smith 	asprintf(&fn, "%d-procs_case/geombc.dat.%d",size,rank+1);
226*16223cb9SCameron Smith 	openfile(fn, "read", &igeombc);
227*16223cb9SCameron Smith 
228*16223cb9SCameron Smith 	std::vector<std::string> headers;
229*16223cb9SCameron Smith 	Gather_Headers(&igeombc, headers);
230*16223cb9SCameron Smith 
231*16223cb9SCameron Smith 	readheader(&igeombc, "size of ilwork array", (void*) &nilwork, &ione,
232*16223cb9SCameron Smith 			"integer", iformat);
233*16223cb9SCameron Smith 	assert(nilwork > 0);
234*16223cb9SCameron Smith 	readheader(&igeombc, "ilwork", (void*) &tmp, &ione,
235*16223cb9SCameron Smith 			"integer", iformat);
236*16223cb9SCameron Smith 	assert(tmp == nilwork);
237*16223cb9SCameron Smith 
238*16223cb9SCameron Smith 	ilwork = (int*) malloc(sizeof(int)*nilwork);
239*16223cb9SCameron Smith 
240*16223cb9SCameron Smith 	readdatablock(&igeombc, "ilwork", ilwork, &nilwork, "integer", iformat);
241*16223cb9SCameron Smith 	closefile(&igeombc, "read");
242*16223cb9SCameron Smith }
243