xref: /phasta/phSolver/common/phasta.cc (revision fe88b52d35537d3736cc3d33d9db7e51353aba7d)
1 //MR CHANGE
2 #define OMPI_SKIP_MPICXX 1
3 //MR CHANGE END
4 #include <mpi.h>
5 #include <stdio.h>
6 #include <stdlib.h>
7 #include <string.h>
8 
9 #include <sys/types.h>
10 #include <sys/stat.h>
11 
12 #include "common_c.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 //MR CHANGE
40 extern void setIOparam(); // For SyncIO
41 //MR CHANGE END
42 
43 int myrank; /* made file global for ease in debugging */
44 
45 void
46 catchDebugger() {
47     while (1) {
48       int debuggerPresent=0;
49       int fakeSTOP = 1; // please stop HERE and assign as next line
50       // assign or set debuggerPresent=1
51       if(debuggerPresent) {
52         break;
53       }
54     }
55 }
56 
57 // some useful debugging functions
58 
59 void
60 pdarray( void* darray , int start, int end ) {
61     for( int i=start; i < end; i++ ){
62         cout << ((double*)darray)[i] << endl;
63     }
64 }
65 
66 void
67 piarray( void* iarray , int start, int end ) {
68     for( int i=start; i < end; i++ ){
69         cout << ((int*)iarray)[i] << endl;
70     }
71 }
72 
73 extern "C" int
74 phasta( int argc,
75         char *argv[] ) {
76 
77     int size,ierr;
78     char inpfilename[100];
79     char* pauseDebugger = getenv("catchDebugger");
80     //cout << "pauseDebugger" << pauseDebugger << endl;
81     int initialized;
82     MPI_Initialized(&initialized);
83     if( !initialized )
84       MPI_Init(&argc,&argv);
85     MPI_Comm_size (MPI_COMM_WORLD, &size);
86     MPI_Comm_rank (MPI_COMM_WORLD, &myrank);
87 
88     workfc.numpe = size;
89     workfc.myrank = myrank;
90 
91 #if (defined WIN32)
92     if(argc > 2 ){
93       catchDebugger();
94     }
95 #endif
96 #if (1) // ALWAYS ( defined LAUNCH_GDB ) && !( defined WIN32 )
97 
98     if ( pauseDebugger ) {
99 
100         int parent_pid = getpid();
101         int gdb_child = fork();
102         cout << "gdb_child" << gdb_child << endl;
103 
104         if( gdb_child == 0 ) {
105 
106             cout << "Debugger Process initiating" << endl;
107             stringstream exec_string;
108 
109 #if ( defined decalp )
110             exec_string <<"xterm -e idb "
111                         << " -pid "<< parent_pid <<" "<< argv[0] << endl;
112 #endif
113 #if ( defined LINUX )
114             exec_string <<"xterm -e gdb"
115                         << " -pid "<< parent_pid <<" "<< argv[0] << endl;
116 #endif
117 #if ( defined SUN4 )
118             exec_string <<"xterm -e dbx "
119                         << " - "<< parent_pid <<" "<< argv[0] << endl;
120 #endif
121 #if ( defined IRIX )
122             exec_string <<"xterm -e dbx "
123                         << " -p "<< parent_pid <<" "<< argv[0] << endl;
124 #endif
125             string s = exec_string.str();
126             system( s.c_str() );
127             exit(0);
128         }
129         catchDebugger();
130     }
131 
132 #endif
133 
134     /* Input data  */
135     if(argc > 1 ){
136         strcpy(inpfilename,argv[1]);
137     } else {
138         strcpy(inpfilename,"solver.inp");
139     }
140     ierr = input_fform(inpfilename);
141     if(!ierr){
142       sprintf(inpfilename,"%d-procs_case/",size);
143       if( chdir( inpfilename ) ) {
144         cerr << "could not change to the problem directory "
145           << inpfilename << endl;
146         return 1;
147       }
148       MPI_Barrier(MPI_COMM_WORLD);
149 
150 //MR CHANGE
151         setIOparam();
152 //MR CHANGE END
153 
154         input();
155         /* now we can start the solver */
156         proces();
157     }
158     else{
159         printf("error during reading ascii input \n");
160     }
161 
162 //MR CHANGE
163 
164     MPI_Barrier(MPI_COMM_WORLD);
165     if ( myrank == 0 ) {
166       printf("phasta.cc - last call before finalize!\n");
167     }
168 //MR CHANGE
169 
170     MPI_Finalize();
171     return 0;
172 }
173