ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-2.0/src/io/basic_ifstrstream.hpp
(Generate patch)

Comparing trunk/OOPSE-2.0/src/io/basic_ifstrstream.hpp (file contents):
Revision 1570 by tim, Fri Oct 15 07:52:47 2004 UTC vs.
Revision 1583 by tim, Fri Oct 15 21:10:50 2004 UTC

# Line 33 | Line 33
33   #ifndef IO_IFSTRSTREAM_HPP
34   #define IO_IFSTRSTREAM_HPP
35  
36 + #include <cassert>
37 + #include <cstring>
38   #include <fstream>
39   #include <sstream>
40  
41 + #ifdef IS_MPI
42 + #include <mpi.h>
43 + #endif
44 +
45   namespace oopse {
46   using namespace std;
47   /**
# Line 50 | Line 56 | using namespace std;
56   *       char buffer[MAXLEN];
57   *       ifstrstream in;
58   *       in.open("Shapes.frc");
59 < *       if (in.is_open) {
59 > *       if (in.is_open()) {
60   *           in.getline(buffer, MAXLEN);
61   *       }
62   *       in.close();
# Line 102 | Line 108 | class basic_ifstrstream : public basic_istream<_CharT,
108          /**
109           * Explicit constructor
110           * @filename String containing the name of the file to be opened
111 <         * @mode Flags describing the requested i/o mode for the file, default value is ios_base::in        
111 >         * @mode Flags describing the requested i/o mode for the file, default value is ios_base::in      
112 >         * @checkFilename Flags indicating checking the file name in parallel
113           */
114 <        explicit basic_ifstrstream(const char* filename, ios_base::openmode mode = ios_base::in)
114 >        explicit basic_ifstrstream(const char* filename, ios_base::openmode mode = ios_base::in, bool checkFilename = false)
115              : basic_ios<_CharT, _Traits>(),  basic_istream<_CharT, _Traits>(0),
116                internalBuf_(NULL), isRead(false) {
117  
118 <           isRead =  internalOpen(filename,  mode | ios_base::in);
118 >           isRead =  internalOpen(filename,  mode | ios_base::in, checkFilename);
119           }
120  
121          /**
122 <         *
122 >         * virtual destructor will close the file(in single mode) and clear the stream buffer
123           */
124          ~basic_ifstrstream(){
125              close();
# Line 126 | Line 133 | class basic_ifstrstream : public basic_istream<_CharT,
133           * brocasting, every nodes fall back to stringstream (parallel mode).
134           * @filename String containing the name of the file to be opened
135           * @mode Flags describing the requested i/o mode for the file
136 +         * @checkFilename Flags indicating checking the file name in parallel
137           */
138 <        void open(const char* filename, ios_base::openmode mode = ios_base::in){
138 >        void open(const char* filename, ios_base::openmode mode = ios_base::in, bool checkFilename = false){
139  
140              if (!isRead ) {
141 <                isRead = internalOpen();
141 >                isRead = internalOpen(filename, mode, checkFilename);
142              }
143          }
144  
# Line 168 | Line 176 | class basic_ifstrstream : public basic_istream<_CharT,
176           * parallel mode) associated with the stream.
177           */
178          _Buf* rdbuf() const{
179 <            return const_cast<_Buf*>(internalBuf_);
179 >            return static_cast<_Buf*>(internalBuf_);
180          }
181  
182      private:
# Line 181 | Line 189 | class basic_ifstrstream : public basic_istream<_CharT,
189           * @mode Flags describing the requested i/o mode for the file
190           * @todo use try - catch syntax to make the program more readable
191           */
192 <        bool internalOpen(const char* filename, ios_base::openmode mode){
192 >        bool internalOpen(const char* filename, ios_base::openmode mode, bool checkFilename){
193  
194   #ifdef IS_MPI        
195              int commStatus;
# Line 190 | Line 198 | class basic_ifstrstream : public basic_istream<_CharT,
198              int filenameLen;
199              int diffFilename;
200              int error;
201 +            int myRank;
202 +            int masterNode;
203  
204 <            if (worldRank == masterNode) {
204 >            commStatus = MPI_Comm_rank(MPI_COMM_WORLD, &myRank);
205 >            masterNode = 0;
206 >            
207 >            if (myRank == masterNode) {
208  
209 <                //check the filename is the same
197 <                filenameLen = strlen(filename);
198 <                commStatus = MPI_Bcast(filenameLen, 1, MPI_INT, masterNode, MPI_COMM_WORLD);    
199 <                commStatus = MPI_Bcast(filename, filenameLen, MPI_CHAR, masterNode, MPI_COMM_WORLD);    
209 >                if (checkFilename) {
210  
211 <                diffFilename = 0;
212 <                commStatus = MPI_Allreduce(sameFilename, error, 1,  MPI_INT,  MPI_COMM_WORLD);            
211 >                    //check the filename is the same
212 >                    filenameLen = strlen(filename);
213 >                    commStatus = MPI_Bcast(&filenameLen, 1, MPI_INT, masterNode, MPI_COMM_WORLD);    
214 >                    commStatus = MPI_Bcast((void*)filename, filenameLen, MPI_CHAR, masterNode, MPI_COMM_WORLD);    
215  
216 <                //if file names are different just return false
217 <                if (error > 0)
216 >                    diffFilename = 0;
217 >                    commStatus = MPI_Allreduce(&diffFilename, &error, 1,  MPI_INT, MPI_SUM,  MPI_COMM_WORLD);            
218 >
219 >                    //if file names are different just return false
220 >                    if (error > 0)
221                          return false;  
222 +                }
223                  
224 <                ifstream fin(filename, mod);
224 >                ifstream fin(filename, mode);
225                  basic_stringbuf<_CharT, _Traits, _Alloc>* sbuf;
226  
227                  if (fin.is_open()) {
228                      
229 <                    fin.in.seekg(0, ios::end);
229 >                    fin.seekg(0, ios::end);
230                      fileSize = fin.tellg();
231 <
231 >                    fin.seekg(0, ios::beg);
232 >                    
233                      // '\0' need one more char
234                      fbuf = new char[fileSize+1];
235                      
# Line 256 | Line 273 | class basic_ifstrstream : public basic_istream<_CharT,
273                
274               } else{
275                      //check file name
276 <                    commStatus = MPI_Bcast(filenameLen, 1, MPI_INT, masterNode, MPI_COMM_WORLD);    
276 >                    if (checkFilename) {
277 >                        commStatus = MPI_Bcast(&filenameLen, 1, MPI_INT, masterNode, MPI_COMM_WORLD);    
278  
279 <                    char * masterFilename = new char[filenameLen];
280 <                    commStatus = MPI_Bcast(masterFilename, filenameLen, MPI_CHAR, masterNode, MPI_COMM_WORLD);    
281 <    
282 <                    if( strcmp(masterFilename, filename) == 0)
283 <                        diffFilename = 0;
284 <                    else
285 <                        diffFilename = 1;
279 >                        char * masterFilename = new char[filenameLen];
280 >                        commStatus = MPI_Bcast(masterFilename, filenameLen, MPI_CHAR, masterNode, MPI_COMM_WORLD);    
281 >        
282 >                        if( strcmp(masterFilename, filename) == 0)
283 >                            diffFilename = 0;
284 >                        else
285 >                            diffFilename = 1;
286  
287 <                    delete masterFilename;
288 <                    
289 <                    commStatus = MPI_Allreduce(sameFilename, error, 1,  MPI_INT,  MPI_COMM_WORLD);    
287 >                        delete masterFilename;
288 >                        
289 >                        commStatus = MPI_Allreduce(&diffFilename, &error, 1,  MPI_INT,  MPI_SUM, MPI_COMM_WORLD);    
290  
291 <                    if (error > 0)
292 <                        return false;                        
293 <
291 >                        if (error > 0)
292 >                            return false;                        
293 >                    }
294                      //get file size
295                      commStatus = MPI_Bcast(&fileSize, 1, MPI_LONG, masterNode, MPI_COMM_WORLD);  
296  
# Line 322 | Line 340 | typedef basic_istringstream<char, char_traits<char>, a
340          bool isRead;                                                                    /** file opened flag */
341   };
342  
343 < typedef basic_istringstream<char, char_traits<char>, allocator<char> > ifstringstream;
343 > typedef basic_ifstrstream<char, char_traits<char>, allocator<char> > ifstrstream;
344   }//namespace oopse
345   #endif //IO_IFSTRSTREAM_HPP

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines