xref: /phasta/checkphasta/checkphasta.cpp (revision dc9538425e8e2132ad0d36a4b7a0472f4a4e7110)
116223cb9SCameron Smith #include <cstdio>
216223cb9SCameron Smith #include <cstdlib>
316223cb9SCameron Smith #include <cassert>
416223cb9SCameron Smith #include <cstring>
516223cb9SCameron Smith #include <cmath>
616223cb9SCameron Smith #include <string>
716223cb9SCameron Smith #include <vector>
816223cb9SCameron Smith #include <algorithm>
916223cb9SCameron Smith #include <set>
1016223cb9SCameron Smith 
1116223cb9SCameron Smith #include <sys/types.h>
1216223cb9SCameron Smith #include <sys/stat.h>
1316223cb9SCameron Smith #include <dirent.h>
1416223cb9SCameron Smith 
1516223cb9SCameron Smith #define OMPI_SKIP_MPICXX 1
1616223cb9SCameron Smith #include <mpi.h>
1716223cb9SCameron Smith 
18bafa09feSCameron Smith #include "phIO.h"
1916223cb9SCameron Smith 
2016223cb9SCameron Smith char read_solution(double** solutiono, int* size, int* nshgo, int* ndofo,
21bafa09feSCameron Smith 		int nump, int rank, int timestep, int nSyncFiles, char* casedir);
22bafa09feSCameron Smith std::set<int>* find_timesteps(char* casedir, int nSyncFiles);
23bafa09feSCameron Smith double compare_solution(char* lpath, char* rpath,
24bafa09feSCameron Smith     int timestep, int nump, int nSyncFiles);
25bafa09feSCameron Smith const char* getRestartName(int nSyncFiles);
2616223cb9SCameron Smith 
2716223cb9SCameron Smith int main(int argc, char** argv)
2816223cb9SCameron Smith {
2916223cb9SCameron Smith        	int rank;
3016223cb9SCameron Smith 	int size;
3116223cb9SCameron Smith 	MPI_Init(&argc, &argv);
3216223cb9SCameron Smith 
3316223cb9SCameron Smith 	MPI_Comm_rank(MPI_COMM_WORLD, &rank);
3416223cb9SCameron Smith 	MPI_Comm_size(MPI_COMM_WORLD, &size);
3516223cb9SCameron Smith 
36c1ff69edSCameron Smith 	if(argc != 5) {
37bafa09feSCameron Smith           fprintf(stderr, "argc %d\n", argc);
38bafa09feSCameron Smith           fprintf(stderr,
39c1ff69edSCameron Smith               "Usage: %s <left> <right> <numSyncFiles> <tolerance>\n"
40bafa09feSCameron Smith               "where <left> and <right> are different"
41bafa09feSCameron Smith               "N-procs_case directories\n", argv[0]);
42bafa09feSCameron Smith           MPI_Finalize();
43bafa09feSCameron Smith           return 1;
44bafa09feSCameron Smith         }
4516223cb9SCameron Smith 	char* lpath = argv[1];
4616223cb9SCameron Smith 	char* rpath = argv[2];
47bafa09feSCameron Smith         int nSyncFiles = atoi(argv[3]);
48c1ff69edSCameron Smith         double tolerance = atof(argv[4]);
4916223cb9SCameron Smith 
5016223cb9SCameron Smith 	int ndof;
5116223cb9SCameron Smith 	int nshg;
5216223cb9SCameron Smith 	int solsize;
5316223cb9SCameron Smith 	double* solution;
5416223cb9SCameron Smith 
55bafa09feSCameron Smith 	std::set<int>* l_timesteps = find_timesteps(lpath, nSyncFiles);
56bafa09feSCameron Smith 	std::set<int>* r_timesteps = find_timesteps(rpath, nSyncFiles);
5716223cb9SCameron Smith 	std::set<int>* timesteps_to_check = new std::set<int>;
5816223cb9SCameron Smith 	std::set_intersection(l_timesteps->begin(), l_timesteps->end(),
5916223cb9SCameron Smith 			r_timesteps->begin(), r_timesteps->end(),
6016223cb9SCameron Smith 			std::inserter(*timesteps_to_check, timesteps_to_check->begin()));
6116223cb9SCameron Smith         delete l_timesteps;
6216223cb9SCameron Smith         delete r_timesteps;
6316223cb9SCameron Smith 	if(rank == 0)
6416223cb9SCameron Smith 		printf("Found %d common timesteps\n",
6516223cb9SCameron Smith 			       	timesteps_to_check->size());
6616223cb9SCameron Smith #ifdef DBGONLY
67bafa09feSCameron Smith 	read_solution(&solution, &solsize, &nshg, &ndof,
68bafa09feSCameron Smith             size, rank, 0, numSyncFiles, "./");
6916223cb9SCameron Smith 	printf("nshg: %d, ndof: %d\n", nshg, ndof);
7016223cb9SCameron Smith 	assert(solsize == ndof*nshg);
7116223cb9SCameron Smith #endif
7216223cb9SCameron Smith 	double maxerror = 0.0;
7316223cb9SCameron Smith 	double error;
7416223cb9SCameron Smith 	double gblmaxerror;
7516223cb9SCameron Smith 	for(std::set<int>::iterator i = timesteps_to_check->begin();
7616223cb9SCameron Smith 			i!=timesteps_to_check->end();i++)
7716223cb9SCameron Smith 	{
78bafa09feSCameron Smith 		error = compare_solution(lpath, rpath, *i, size, nSyncFiles);
7916223cb9SCameron Smith 		if(error>maxerror) maxerror = error;
8016223cb9SCameron Smith 	}
8116223cb9SCameron Smith         delete timesteps_to_check;
8216223cb9SCameron Smith 	MPI_Barrier(MPI_COMM_WORLD);
8316223cb9SCameron Smith 	MPI_Reduce(&maxerror, &gblmaxerror, 1, MPI_DOUBLE, MPI_MAX, 0,
8416223cb9SCameron Smith 		MPI_COMM_WORLD);
8516223cb9SCameron Smith 	if(rank == 0) printf("Maximum difference across all timesteps: %e\n",
8616223cb9SCameron Smith 			gblmaxerror);
8716223cb9SCameron Smith 	MPI_Finalize();
88c1ff69edSCameron Smith         return (gblmaxerror > tolerance);
8916223cb9SCameron Smith }
90bafa09feSCameron Smith double compare_solution(char* lpath, char* rpath, int timestep, int nump, int nSyncFiles)
9116223cb9SCameron Smith {
9216223cb9SCameron Smith 	int rank;
9316223cb9SCameron Smith 	MPI_Comm_rank(MPI_COMM_WORLD, &rank);
9416223cb9SCameron Smith 	double* lsol;
9516223cb9SCameron Smith 	double* rsol;
9616223cb9SCameron Smith 	int lsize;
9716223cb9SCameron Smith 	int rsize;
9816223cb9SCameron Smith 
99bafa09feSCameron Smith 	read_solution(&lsol, &lsize, NULL, NULL,
100bafa09feSCameron Smith             nump, rank, timestep, nSyncFiles, lpath);
101bafa09feSCameron Smith 	read_solution(&rsol, &rsize, NULL, NULL,
102bafa09feSCameron Smith             nump, rank, timestep, nSyncFiles, rpath);
10316223cb9SCameron Smith 
10416223cb9SCameron Smith 	double maxdiff=0.0;
10516223cb9SCameron Smith 	double gblmaxdiff;
10616223cb9SCameron Smith 	if(lsize != rsize)
10716223cb9SCameron Smith 	{
10816223cb9SCameron Smith 		printf("Error: Solution sizes different: %d, %d\n",
10916223cb9SCameron Smith 				lsize, rsize);
11016223cb9SCameron Smith 		assert(lsize == rsize);
11116223cb9SCameron Smith 	}
11216223cb9SCameron Smith 	for(int i=0;i<lsize;i++)
11316223cb9SCameron Smith 	{
11416223cb9SCameron Smith 		double diff = fabs(lsol[i]-rsol[i]);
11516223cb9SCameron Smith 		if(diff > maxdiff) maxdiff = diff;
11616223cb9SCameron Smith 	}
11716223cb9SCameron Smith         free(lsol);
11816223cb9SCameron Smith         free(rsol);
11916223cb9SCameron Smith 	MPI_Reduce(&maxdiff, &gblmaxdiff, 1, MPI_DOUBLE, MPI_MAX, 0,
12016223cb9SCameron Smith 		       	MPI_COMM_WORLD);
12116223cb9SCameron Smith 	MPI_Barrier(MPI_COMM_WORLD); //TODO: debugging only
12216223cb9SCameron Smith 	if(rank == 0)
12316223cb9SCameron Smith 	printf("Timestep: %d, maximum difference: %e\n", timestep, gblmaxdiff);
12416223cb9SCameron Smith 	return gblmaxdiff;
12516223cb9SCameron Smith 
12616223cb9SCameron Smith }
12716223cb9SCameron Smith char read_solution(double** solutiono, int* size, int* nshgo, int* ndofo,
128bafa09feSCameron Smith 		int nump, int rank, int timestep, int nSyncFiles, char* casedir)
12916223cb9SCameron Smith {
13016223cb9SCameron Smith 	int iarray[10];
13116223cb9SCameron Smith         const char* iformat = "binary";
13216223cb9SCameron Smith         int ithree=3;
13316223cb9SCameron Smith         int igeombc;
13416223cb9SCameron Smith         char* fn;
13516223cb9SCameron Smith 	int nshg;
13616223cb9SCameron Smith 	int ndof;
13716223cb9SCameron Smith 	double* solution;
138bafa09feSCameron Smith         const char* rname = getRestartName(nSyncFiles);
1390acf3134SCameron Smith         asprintf(&fn,"%s/%s.%d.",casedir,rname,timestep);
140bafa09feSCameron Smith         phio_fp fp;
141bafa09feSCameron Smith         phio_openfile_read(fn,&nSyncFiles,&fp);
142bafa09feSCameron Smith 	phio_readheader(fp, "solution", (void*) iarray, &ithree, "integer", iformat);
14316223cb9SCameron Smith 	nshg = iarray[0];
14416223cb9SCameron Smith 	ndof = iarray[1];
14516223cb9SCameron Smith 	if(size != NULL)
14616223cb9SCameron Smith 		*size = nshg*ndof;
14716223cb9SCameron Smith 	solution = (double*) malloc(sizeof(double)*nshg*ndof);
148bafa09feSCameron Smith 	phio_readdatablock(fp, "solution", solution, size, "double", iformat);
149bafa09feSCameron Smith 	phio_closefile_read(fp);
15016223cb9SCameron Smith 	if(solutiono != NULL)
15116223cb9SCameron Smith 		*solutiono = solution;
15216223cb9SCameron Smith 	if(nshgo != NULL)
15316223cb9SCameron Smith 		*nshgo = nshg;
15416223cb9SCameron Smith 	if(ndofo != NULL)
15516223cb9SCameron Smith 		*ndofo = ndof;
15616223cb9SCameron Smith 	free(fn);
15716223cb9SCameron Smith 	return(0);
15816223cb9SCameron Smith }
15916223cb9SCameron Smith 
160bafa09feSCameron Smith std::set<int>* find_timesteps(char* casedir, int nSyncFiles)
16116223cb9SCameron Smith {
16216223cb9SCameron Smith 	char* path;
16316223cb9SCameron Smith 	DIR* dir;
16416223cb9SCameron Smith 	struct dirent* d;
16516223cb9SCameron Smith 	int part, ts;
16616223cb9SCameron Smith 	std::set<int>* step_list = new std::set<int>;
16716223cb9SCameron Smith 
16816223cb9SCameron Smith         asprintf(&path, "%s", casedir);
16916223cb9SCameron Smith 	dir = opendir(path);
17016223cb9SCameron Smith 	if(!dir)
17116223cb9SCameron Smith 	{
17216223cb9SCameron Smith 		perror("Error opening case: ");
17316223cb9SCameron Smith 		MPI_Abort(MPI_COMM_WORLD,1);
17416223cb9SCameron Smith 	}
175bafa09feSCameron Smith         const char* rname = getRestartName(nSyncFiles);
176bafa09feSCameron Smith         char* fmt;
177bafa09feSCameron Smith         asprintf(&fmt, "%s.%%d.%%d", rname);
17816223cb9SCameron Smith 	while((d=readdir(dir)))
17916223cb9SCameron Smith 	{
1800acf3134SCameron Smith 		if(sscanf(d->d_name, fmt, &ts, &part)==2)
18116223cb9SCameron Smith 		{
18216223cb9SCameron Smith 			step_list->insert(ts);
18316223cb9SCameron Smith 		}
18416223cb9SCameron Smith 	}
185bafa09feSCameron Smith         free(fmt);
18616223cb9SCameron Smith 	free(path);
187*dc953842SCameron Smith         closedir(dir);
188bafa09feSCameron Smith 	return(step_list);
189bafa09feSCameron Smith }
190bafa09feSCameron Smith 
191bafa09feSCameron Smith const char* getRestartName(int nSyncFiles) {
192bafa09feSCameron Smith   if(0 == nSyncFiles)
193bafa09feSCameron Smith     return "restart";
194bafa09feSCameron Smith   else if(nSyncFiles > 0)
195bafa09feSCameron Smith     return "restart-dat";
196bafa09feSCameron Smith   else {
197bafa09feSCameron Smith     fprintf(stderr,
198bafa09feSCameron Smith         "ERROR: the number of sync-io files must be"
199bafa09feSCameron Smith         "greater than or equal to zero\n");
200bafa09feSCameron Smith     MPI_Abort(MPI_COMM_WORLD, 1);
201bafa09feSCameron Smith     return NULL;
202bafa09feSCameron Smith   }
20316223cb9SCameron Smith }
204