VTK
vtkMPIController.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkMPIController.h
5 
6  Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
7  All rights reserved.
8  See Copyright.txt or http://www.kitware.com/Copyright.htm for details.
9 
10  This software is distributed WITHOUT ANY WARRANTY; without even
11  the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
12  PURPOSE. See the above copyright notice for more information.
13 
14 =========================================================================*/
48 #ifndef vtkMPIController_h
49 #define vtkMPIController_h
50 
51 #include "vtkParallelMPIModule.h" // For export macro
53 // Do not remove this header file. This class contains methods
54 // which take arguments defined in vtkMPICommunicator.h by
55 // reference.
56 #include "vtkMPICommunicator.h" // Needed for direct access to communicator
57 
58 class vtkIntArray;
59 
61 {
62 
63 public:
64 
65  static vtkMPIController *New();
67  void PrintSelf(ostream& os, vtkIndent indent);
68 
70 
78  virtual void Initialize(int* argc, char*** argv)
79  { this->Initialize(argc, argv, 0); }
81 
82  virtual void Initialize(int* vtkNotUsed(argc), char*** vtkNotUsed(argv),
83  int initializedExternally);
84 
87  virtual void Initialize();
88 
91  virtual void Finalize() { this->Finalize(0); }
92 
93  virtual void Finalize(int finalizedExternally);
94 
97  virtual void SingleMethodExecute();
98 
102  virtual void MultipleMethodExecute();
103 
106  virtual void CreateOutputWindow();
107 
110  static char* ErrorString(int err);
111 
112 
118  void SetCommunicator(vtkMPICommunicator* comm);
119 
121 
122  virtual vtkMPIController *PartitionController(int localColor, int localKey);
123 
124 //BTX
125 
127 
133  int NoBlockSend(const int* data, int length, int remoteProcessId, int tag,
135  { return ((vtkMPICommunicator*)this->Communicator)->NoBlockSend
136  (data ,length, remoteProcessId, tag, req); }
137  int NoBlockSend(const unsigned long* data, int length, int remoteProcessId,
138  int tag, vtkMPICommunicator::Request& req)
139  { return ((vtkMPICommunicator*)this->Communicator)->NoBlockSend
140  (data, length, remoteProcessId, tag, req); }
141  int NoBlockSend(const char* data, int length, int remoteProcessId,
142  int tag, vtkMPICommunicator::Request& req)
143  { return ((vtkMPICommunicator*)this->Communicator)->NoBlockSend
144  (data, length, remoteProcessId, tag, req); }
145  int NoBlockSend( const unsigned char* data, int length, int remoteProcessId,
146  int tag, vtkMPICommunicator::Request& req )
147  { return ((vtkMPICommunicator*)this->Communicator)->NoBlockSend
148  (data, length, remoteProcessId, tag, req);}
149  int NoBlockSend(const float* data, int length, int remoteProcessId,
150  int tag, vtkMPICommunicator::Request& req)
151  { return ((vtkMPICommunicator*)this->Communicator)->NoBlockSend
152  (data, length, remoteProcessId, tag, req); }
153  int NoBlockSend(const double* data, int length, int remoteProcessId,
154  int tag, vtkMPICommunicator::Request& req)
155  { return ((vtkMPICommunicator*)this->Communicator)->NoBlockSend
156  (data, length, remoteProcessId, tag, req); }
157 #ifdef VTK_USE_64BIT_IDS
158  int NoBlockSend(const vtkIdType* data, int length, int remoteProcessId,
159  int tag, vtkMPICommunicator::Request& req)
160  { return ((vtkMPICommunicator*)this->Communicator)->NoBlockSend
161  (data, length, remoteProcessId, tag, req); }
162 #endif
163 
164 
166 
171  int NoBlockReceive(int* data, int length, int remoteProcessId,
172  int tag, vtkMPICommunicator::Request& req)
173  { return ((vtkMPICommunicator*)this->Communicator)->NoBlockReceive
174  (data, length, remoteProcessId, tag, req); }
175  int NoBlockReceive(unsigned long* data, int length,
176  int remoteProcessId, int tag,
178  { return ((vtkMPICommunicator*)this->Communicator)->NoBlockReceive
179  (data, length, remoteProcessId, tag, req); }
180  int NoBlockReceive(char* data, int length, int remoteProcessId,
181  int tag, vtkMPICommunicator::Request& req)
182  { return ((vtkMPICommunicator*)this->Communicator)->NoBlockReceive
183  (data, length, remoteProcessId, tag, req); }
184  int NoBlockReceive(unsigned char* data, int length, int remoteProcessId,
185  int tag, vtkMPICommunicator::Request& req)
186  { return ((vtkMPICommunicator*)this->Communicator)->NoBlockReceive
187  (data, length, remoteProcessId, tag, req); }
188  int NoBlockReceive(float* data, int length, int remoteProcessId,
189  int tag, vtkMPICommunicator::Request& req)
190  { return ((vtkMPICommunicator*)this->Communicator)->NoBlockReceive
191  (data, length, remoteProcessId, tag, req); }
192  int NoBlockReceive(double* data, int length, int remoteProcessId,
193  int tag, vtkMPICommunicator::Request& req)
194  { return ((vtkMPICommunicator*)this->Communicator)->NoBlockReceive
195  (data, length, remoteProcessId, tag, req); }
196 #ifdef VTK_USE_64BIT_IDS
197  int NoBlockReceive(vtkIdType* data, int length, int remoteProcessId,
198  int tag, vtkMPICommunicator::Request& req)
199  { return ((vtkMPICommunicator*)this->Communicator)->NoBlockReceive
200  (data, length, remoteProcessId, tag, req); }
201 #endif
202 
203 
205 
212  int Iprobe(int source, int tag, int* flag, int* actualSource)
213  { return ((vtkMPICommunicator*)this->Communicator)->Iprobe(
214  source, tag, flag, actualSource); }
215  int Iprobe(int source, int tag, int* flag, int* actualSource,
216  int* type, int* size)
217  { return ((vtkMPICommunicator*)this->Communicator)->Iprobe(
218  source, tag, flag, actualSource, type, size); }
219  int Iprobe(int source, int tag, int* flag, int* actualSource,
220  unsigned long* type, int* size)
221  { return ((vtkMPICommunicator*)this->Communicator)->Iprobe(
222  source, tag, flag, actualSource, type, size); }
223  int Iprobe(int source, int tag, int* flag, int* actualSource,
224  const char* type, int* size)
225  { return ((vtkMPICommunicator*)this->Communicator)->Iprobe(
226  source, tag, flag, actualSource, type, size); }
227  int Iprobe(int source, int tag, int* flag, int* actualSource,
228  float* type, int* size)
229  { return ((vtkMPICommunicator*)this->Communicator)->Iprobe(
230  source, tag, flag, actualSource, type, size); }
231  int Iprobe(int source, int tag, int* flag, int* actualSource,
232  double* type, int* size)
233  { return ((vtkMPICommunicator*)this->Communicator)->Iprobe(
234  source, tag, flag, actualSource, type, size); }
236 
238 
241  int WaitAll(const int count, vtkMPICommunicator::Request requests[])
242  { return ((vtkMPICommunicator*)this->Communicator)->WaitAll(count,requests);}
244 
246 
250  int WaitAny(const int count,vtkMPICommunicator::Request requests[], int& idx)
251  {return ((vtkMPICommunicator*)this->Communicator)->WaitAny(count,requests,idx);}
253 
255 
258  int WaitSome(
259  const int count, vtkMPICommunicator::Request requests[],
260  vtkIntArray *completed );
262 
265  bool TestAll(const int count, vtkMPICommunicator::Request requests[] );
266 
271  bool TestAny(const int count,vtkMPICommunicator::Request requests[],int &idx);
272 
274 
277  bool TestSome(const int count,vtkMPICommunicator::Request requests[],
278  vtkIntArray *completed );
279 //ETX
280  static const char* GetProcessorName();
282 
284 
286  static void SetUseSsendForRMI(int use_send)
287  { vtkMPIController::UseSsendForRMI = (use_send != 0)? 1: 0; }
289 //BTX
290 protected:
292  ~vtkMPIController();
294 
295  // Set the communicator to comm and call InitializeNumberOfProcesses()
296  void InitializeCommunicator(vtkMPICommunicator* comm);
297 
298  // Duplicate the current communicator, creating RMICommunicator
299  void InitializeRMICommunicator();
300 
302 
305  virtual void TriggerRMIInternal(int remoteProcessId,
306  void* arg, int argLength, int rmiTag, bool propagate);
308 
309  // MPI communicator created when Initialize() called.
310  // This is a copy of MPI_COMM_WORLD but uses a new
311  // context, i.e. even if the tags are the same, the
312  // RMI messages will not interfere with user level
313  // messages.
315 
316  friend class vtkMPIOutputWindow;
317 
318  // Initialize only once.
319  static int Initialized;
320 
321  static char ProcessorName[];
322 
324 
325  static int UseSsendForRMI;
326 private:
327  vtkMPIController(const vtkMPIController&); // Not implemented.
328  void operator=(const vtkMPIController&); // Not implemented.
329 //ETX
330 };
332 
333 
334 #endif
335 
336 
int NoBlockSend(const double *data, int length, int remoteProcessId, int tag, vtkMPICommunicator::Request &req)
int NoBlockSend(const unsigned char *data, int length, int remoteProcessId, int tag, vtkMPICommunicator::Request &req)
int NoBlockSend(const unsigned long *data, int length, int remoteProcessId, int tag, vtkMPICommunicator::Request &req)
int WaitAll(const int count, vtkMPICommunicator::Request requests[])
int WaitAny(const int count, vtkMPICommunicator::Request requests[], int &idx)
virtual void Finalize()=0
Class for creating user defined MPI communicators.
int Iprobe(int source, int tag, int *flag, int *actualSource, const char *type, int *size)
virtual void TriggerRMIInternal(int remoteProcessId, void *arg, int argLength, int rmiTag, bool propagate)
int vtkIdType
Definition: vtkType.h:247
virtual void MultipleMethodExecute()=0
static int UseSsendForRMI
int NoBlockReceive(unsigned char *data, int length, int remoteProcessId, int tag, vtkMPICommunicator::Request &req)
static void SetUseSsendForRMI(int use_send)
dynamic, self-adjusting array of int
Definition: vtkIntArray.h:49
int NoBlockSend(const int *data, int length, int remoteProcessId, int tag, vtkMPICommunicator::Request &req)
virtual void Initialize(int *vtkNotUsed(argc), char ***vtkNotUsed(argv))=0
int NoBlockReceive(unsigned long *data, int length, int remoteProcessId, int tag, vtkMPICommunicator::Request &req)
a simple class to control print indentation
Definition: vtkIndent.h:38
static int Initialized
A subgroup of processes from a communicator.
Process communication using MPI.
static vtkMPICommunicator * WorldRMICommunicator
int NoBlockReceive(float *data, int length, int remoteProcessId, int tag, vtkMPICommunicator::Request &req)
virtual void Initialize(int *argc, char ***argv)
int Iprobe(int source, int tag, int *flag, int *actualSource, double *type, int *size)
virtual void CreateOutputWindow()=0
int NoBlockSend(const char *data, int length, int remoteProcessId, int tag, vtkMPICommunicator::Request &req)
boost::graph_traits< vtkGraph * >::vertex_descriptor source(boost::graph_traits< vtkGraph * >::edge_descriptor e, vtkGraph *)
int Iprobe(int source, int tag, int *flag, int *actualSource, unsigned long *type, int *size)
int NoBlockReceive(double *data, int length, int remoteProcessId, int tag, vtkMPICommunicator::Request &req)
virtual void Finalize()
int NoBlockReceive(char *data, int length, int remoteProcessId, int tag, vtkMPICommunicator::Request &req)
void PrintSelf(ostream &os, vtkIndent indent)
int Iprobe(int source, int tag, int *flag, int *actualSource, float *type, int *size)
static int GetUseSsendForRMI()
int NoBlockReceive(int *data, int length, int remoteProcessId, int tag, vtkMPICommunicator::Request &req)
static vtkObject * New()
int NoBlockSend(const float *data, int length, int remoteProcessId, int tag, vtkMPICommunicator::Request &req)
int Iprobe(int source, int tag, int *flag, int *actualSource, int *type, int *size)
virtual void SingleMethodExecute()=0
virtual vtkMultiProcessController * PartitionController(int localColor, int localKey)
virtual vtkMultiProcessController * CreateSubController(vtkProcessGroup *group)
int Iprobe(int source, int tag, int *flag, int *actualSource)
Multiprocessing communication superclass.
#define VTKPARALLELMPI_EXPORT