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