xref: /phasta/phSolver/common/phasta.cc (revision 595995161822a203c8467e0e4a253d7bd7d6df32)
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 <strstream>
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 void Partition_Problem( int, char[], char[], int );
37 extern "C" void proces();
38 extern "C" void input();
39 extern int input_fform(char inpfname[]);
40 //MR CHANGE
41 extern void setIOparam(); // For SyncIO
42 //MR CHANGE END
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 extern "C" int
75 phasta( int argc,
76         char *argv[] ) {
77 
78     int size,ierr;
79     char inpfilename[100];
80     char* pauseDebugger = getenv("catchDebugger");
81     //cout << "pauseDebugger" << pauseDebugger << endl;
82 
83     MPI_Init(&argc,&argv);
84     MPI_Comm_size (MPI_COMM_WORLD, &size);
85     MPI_Comm_rank (MPI_COMM_WORLD, &myrank);
86 
87     workfc.numpe = size;
88     workfc.myrank = myrank;
89 
90 #if (defined WIN32)
91     if(argc > 2 ){
92       catchDebugger();
93     }
94 #endif
95 #if (1) // ALWAYS ( defined LAUNCH_GDB ) && !( defined WIN32 )
96 
97     if ( pauseDebugger ) {
98 
99         int parent_pid = getpid();
100         int gdb_child = fork();
101         cout << "gdb_child" << gdb_child << endl;
102 
103         if( gdb_child == 0 ) {
104 
105             cout << "Debugger Process initiating" << endl;
106             strstream exec_string;
107 
108 #if ( defined decalp )
109             exec_string <<"xterm -e idb "
110                         << " -pid "<< parent_pid <<" "<< argv[0] << endl;
111 #endif
112 #if ( defined LINUX )
113             exec_string <<"xterm -e gdb"
114                         << " -pid "<< parent_pid <<" "<< argv[0] << endl;
115 #endif
116 #if ( defined SUN4 )
117             exec_string <<"xterm -e dbx "
118                         << " - "<< parent_pid <<" "<< argv[0] << endl;
119 #endif
120 #if ( defined IRIX )
121             exec_string <<"xterm -e dbx "
122                         << " -p "<< parent_pid <<" "<< argv[0] << endl;
123 #endif
124             system( exec_string.str() );
125             exit(0);
126         }
127         catchDebugger();
128     }
129 
130 #endif
131 
132     /* Input data  */
133     if(argc > 1 ){
134         strcpy(inpfilename,argv[1]);
135     } else {
136         strcpy(inpfilename,"solver.inp");
137     }
138     ierr = input_fform(inpfilename);
139     if(!ierr){
140 
141         /* Preprocess data and run the problem  */
142         /* Partition the problem to the correct number of processors */
143         if( size > 1 ) {
144             if( myrank == 0 ) {
145                  Partition_Problem( size, phasta_iotype,
146                                     phasta_iotype, SONFATH );
147                  MPI_Barrier(MPI_COMM_WORLD);
148             } else {
149                  MPI_Barrier(MPI_COMM_WORLD);
150                  sprintf(inpfilename,"%d-procs_case/",size);
151                  if( chdir( inpfilename ) ) {
152                     cerr << "could not change to the problem directory "
153                               << inpfilename << endl;
154                     return 1;
155                  }
156             }
157         }
158 
159 //MR CHANGE
160         setIOparam();
161 //MR CHANGE END
162 
163         input();
164         /* now we can start the solver */
165         proces();
166     }
167     else{
168         printf("error during reading ascii input \n");
169     }
170 
171 //MR CHANGE
172 
173     MPI_Barrier(MPI_COMM_WORLD);
174     if ( myrank == 0 ) {
175       printf("phasta.cc - last call before finalize!\n");
176     }
177 //MR CHANGE
178 
179     MPI_Finalize();
180     return 0;
181 }
182