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