xref: /phasta/phSolver/common/phasta.cc (revision fa773e3ffef2f6750f0d17842585fcd2becb595d)
1 #define OMPI_SKIP_MPICXX 1
2 #include <mpi.h>
3 #include <stdio.h>
4 #include <stdlib.h>
5 #include <string.h>
6 #include <cassert>
7 
8 #include <sys/types.h>
9 #include <sys/stat.h>
10 
11 #include "common_c.h"
12 #include "Input.h"
13 #include "phstream.h"
14 #include "streamio.h"
15 
16 #if !(defined IOSTREAMH)
17 #include <iostream>
18 #include <sstream>
19 using namespace std;
20 #endif
21 
22 #include <FCMangle.h>
23 #define input FortranCInterface_GLOBAL_(input,INPUT)
24 #define proces FortranCInterface_GLOBAL_(proces,PROCES)
25 #define timer FortranCInterface_GLOBAL_(timer,TIMER)
26 
27 #ifdef intel
28 #include <direct.h>
29 #define chdir _chdir
30 #else
31 #include <unistd.h>
32 #endif
33 
34 extern "C" char phasta_iotype[80];
35 char phasta_iotype[80];
36 
37 extern int SONFATH;
38 extern "C" void proces();
39 extern "C" void input();
40 extern int input_fform(phSolver::Input&);
41 extern void setIOparam(); // For SyncIO
42 extern "C" void initPhastaCommonVars();
43 
44 int myrank; /* made file global for ease in debugging */
45 
46 void
47 catchDebugger() {
48     while (1) {
49       int debuggerPresent=0;
50       int fakeSTOP = 1; // please stop HERE and assign as next line
51       // assign or set debuggerPresent=1
52       if(debuggerPresent) {
53         break;
54       }
55     }
56 }
57 
58 // some useful debugging functions
59 
60 void
61 pdarray( void* darray , int start, int end ) {
62     for( int i=start; i < end; i++ ){
63         cout << ((double*)darray)[i] << endl;
64     }
65 }
66 
67 void
68 piarray( void* iarray , int start, int end ) {
69     for( int i=start; i < end; i++ ){
70         cout << ((int*)iarray)[i] << endl;
71     }
72 }
73 
74 namespace {
75   int cdToParent() {
76     if( chdir("..") ) {
77       fprintf(stderr,"could not change to the parent directory\n");
78       return 1;
79     } else {
80       return 0;
81     }
82   }
83   int run(phSolver::Input& ctrl) {
84     int size,ierr;
85     char inpfilename[100];
86     MPI_Comm_size (MPI_COMM_WORLD, &size);
87     MPI_Comm_rank (MPI_COMM_WORLD, &myrank);
88 
89     workfc.numpe = size;
90     workfc.myrank = myrank;
91 
92     initPhastaCommonVars();
93     /* Input data  */
94     ierr = input_fform(ctrl);
95     if(!ierr){
96       sprintf(inpfilename,"%d-procs_case/",size);
97       if( chdir( inpfilename ) ) {
98         cerr << "could not change to the problem directory "
99           << inpfilename << endl;
100         return -1;
101       }
102       MPI_Barrier(MPI_COMM_WORLD);
103       input();
104       /* now we can start the solver */
105       proces();
106     }
107     else{
108       printf("error during reading ascii input \n");
109     }
110     MPI_Barrier(MPI_COMM_WORLD);
111     if ( myrank == 0 ) {
112       printf("phasta.cc - last call before finalize!\n");
113     }
114     if( cdToParent() )
115       return -1;
116     return timdat.lstep;
117   }
118 }
119 
120 int phasta(phSolver::Input& ctrl) {
121   outpar.input_mode = 0; //FIXME magic value for posix
122   outpar.output_mode = 0; //FIXME magic value for posix
123   return run(ctrl);
124 }
125 
126 int phasta(phSolver::Input& ctrl, grstream grs) {
127   assert(grs);
128   outpar.input_mode = -1; //FIXME magic value for streams
129   outpar.output_mode = 1; //FIXME magic value for syncio
130   streamio_set_gr(grs);
131   return run(ctrl);
132 }
133 
134 int phasta(phSolver::Input& ctrl, RStream* rs) {
135   fprintf(stderr, "HEY! if you see this email Cameron and tell him "
136       "to implement %s(...) on line %d of %s "
137       "... returning an error\n", __func__, __LINE__, __FILE__);
138   return -1;
139 }
140 
141 int phasta(phSolver::Input& ctrl, GRStream* grs, RStream* rs) {
142   outpar.input_mode = -1; //FIXME magic value for streams
143   outpar.output_mode = -1; //FIXME magic value for streams
144   assert(grs);
145   assert(rs);
146   streamio_set_gr(grs);
147   streamio_set_r(rs);
148   return run(ctrl);
149 }
150 
151 int phasta( int argc, char *argv[] ) {
152     int size,ierr;
153     char inpfilename[100];
154     char* pauseDebugger = getenv("catchDebugger");
155     MPI_Comm_size (MPI_COMM_WORLD, &size);
156     MPI_Comm_rank (MPI_COMM_WORLD, &myrank);
157 
158     workfc.numpe = size;
159     workfc.myrank = myrank;
160 
161 #if (defined WIN32)
162     if(argc > 2 ){
163       catchDebugger();
164     }
165 #endif
166 #if (1) // ALWAYS ( defined LAUNCH_GDB ) && !( defined WIN32 )
167 
168     if ( pauseDebugger ) {
169 
170         int parent_pid = getpid();
171         int gdb_child = fork();
172         cout << "gdb_child" << gdb_child << endl;
173 
174         if( gdb_child == 0 ) {
175 
176             cout << "Debugger Process initiating" << endl;
177             stringstream exec_string;
178 
179 #if ( defined decalp )
180             exec_string <<"xterm -e idb "
181                         << " -pid "<< parent_pid <<" "<< argv[0] << endl;
182 #endif
183 #if ( defined LINUX )
184             exec_string <<"xterm -e gdb"
185                         << " -pid "<< parent_pid <<" "<< argv[0] << endl;
186 #endif
187 #if ( defined SUN4 )
188             exec_string <<"xterm -e dbx "
189                         << " - "<< parent_pid <<" "<< argv[0] << endl;
190 #endif
191 #if ( defined IRIX )
192             exec_string <<"xterm -e dbx "
193                         << " -p "<< parent_pid <<" "<< argv[0] << endl;
194 #endif
195             string s = exec_string.str();
196             system( s.c_str() );
197             exit(0);
198         }
199         catchDebugger();
200     }
201 
202 #endif
203 
204     /* Input data  */
205     if(argc > 1 ){
206         strcpy(inpfilename,argv[1]);
207     } else {
208         strcpy(inpfilename,"solver.inp");
209     }
210     string defaultConf = ".";
211     const char* path_to_config = getenv("PHASTA_CONFIG");
212     if(path_to_config)
213       defaultConf = path_to_config;
214     defaultConf.append("/input.config");
215     string userConf(inpfilename);
216     phSolver::Input ctrl(userConf, defaultConf);
217     initPhastaCommonVars();
218     ierr = input_fform(ctrl);
219     if(!ierr){
220       sprintf(inpfilename,"%d-procs_case/",size);
221       if( chdir( inpfilename ) ) {
222         cerr << "could not change to the problem directory "
223           << inpfilename << endl;
224         return -1;
225       }
226       MPI_Barrier(MPI_COMM_WORLD);
227       setIOparam();
228       outpar.input_mode = outpar.nsynciofiles; //FIXME this is awful
229       outpar.output_mode = outpar.nsynciofiles; //FIXME this is awful
230       input();
231       /* now we can start the solver */
232       proces();
233     }
234     else{
235         printf("error during reading ascii input \n");
236     }
237     MPI_Barrier(MPI_COMM_WORLD);
238     if ( myrank == 0 ) {
239       printf("phasta.cc - last call before finalize!\n");
240     }
241     if( cdToParent() )
242       return -1;
243     return timdat.lstep;
244 }
245