VTK
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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); }
158 
160 
165  int NoBlockReceive(int* data, int length, int remoteProcessId,
166  int tag, vtkMPICommunicator::Request& req)
167  { return ((vtkMPICommunicator*)this->Communicator)->NoBlockReceive
168  (data, length, remoteProcessId, tag, req); }
169  int NoBlockReceive(unsigned long* data, int length,
170  int remoteProcessId, int tag,
172  { return ((vtkMPICommunicator*)this->Communicator)->NoBlockReceive
173  (data, length, remoteProcessId, tag, req); }
174  int NoBlockReceive(char* data, int length, int remoteProcessId,
175  int tag, vtkMPICommunicator::Request& req)
176  { return ((vtkMPICommunicator*)this->Communicator)->NoBlockReceive
177  (data, length, remoteProcessId, tag, req); }
178  int NoBlockReceive(unsigned char* data, int length, int remoteProcessId,
179  int tag, vtkMPICommunicator::Request& req)
180  { return ((vtkMPICommunicator*)this->Communicator)->NoBlockReceive
181  (data, length, remoteProcessId, tag, req); }
182  int NoBlockReceive(float* data, int length, int remoteProcessId,
183  int tag, vtkMPICommunicator::Request& req)
184  { return ((vtkMPICommunicator*)this->Communicator)->NoBlockReceive
185  (data, length, remoteProcessId, tag, req); }
186  int NoBlockReceive(double* data, int length, int remoteProcessId,
187  int tag, vtkMPICommunicator::Request& req)
188  { return ((vtkMPICommunicator*)this->Communicator)->NoBlockReceive
189  (data, length, remoteProcessId, tag, req); }
190 #ifdef VTK_USE_64BIT_IDS
191  int NoBlockReceive(vtkIdType* data, int length, int remoteProcessId,
192  int tag, vtkMPICommunicator::Request& req)
193  { return ((vtkMPICommunicator*)this->Communicator)->NoBlockReceive
194  (data, length, remoteProcessId, tag, req); }
195 #endif
196 
197 
199 
206  int Iprobe(int source, int tag, int* flag, int* actualSource)
207  { return ((vtkMPICommunicator*)this->Communicator)->Iprobe(
208  source, tag, flag, actualSource); }
209  int Iprobe(int source, int tag, int* flag, int* actualSource,
210  int* type, int* size)
211  { return ((vtkMPICommunicator*)this->Communicator)->Iprobe(
212  source, tag, flag, actualSource, type, size); }
213  int Iprobe(int source, int tag, int* flag, int* actualSource,
214  unsigned long* type, int* size)
215  { return ((vtkMPICommunicator*)this->Communicator)->Iprobe(
216  source, tag, flag, actualSource, type, size); }
217  int Iprobe(int source, int tag, int* flag, int* actualSource,
218  const char* type, int* size)
219  { return ((vtkMPICommunicator*)this->Communicator)->Iprobe(
220  source, tag, flag, actualSource, type, size); }
221  int Iprobe(int source, int tag, int* flag, int* actualSource,
222  float* type, int* size)
223  { return ((vtkMPICommunicator*)this->Communicator)->Iprobe(
224  source, tag, flag, actualSource, type, size); }
225  int Iprobe(int source, int tag, int* flag, int* actualSource,
226  double* type, int* size)
227  { return ((vtkMPICommunicator*)this->Communicator)->Iprobe(
228  source, tag, flag, actualSource, type, size); }
230 
232 
235  int WaitAll(const int count, vtkMPICommunicator::Request requests[])
236  { return ((vtkMPICommunicator*)this->Communicator)->WaitAll(count,requests);}
238 
240 
244  int WaitAny(const int count,vtkMPICommunicator::Request requests[], int& idx)
245  {return ((vtkMPICommunicator*)this->Communicator)->WaitAny(count,requests,idx);}
247 
249 
252  int WaitSome(
253  const int count, vtkMPICommunicator::Request requests[],
254  vtkIntArray *completed );
256 
259  bool TestAll(const int count, vtkMPICommunicator::Request requests[] );
260 
265  bool TestAny(const int count,vtkMPICommunicator::Request requests[],int &idx);
266 
268 
271  bool TestSome(const int count,vtkMPICommunicator::Request requests[],
272  vtkIntArray *completed );
273 //ETX
274  static const char* GetProcessorName();
276 
278 
280  static void SetUseSsendForRMI(int use_send)
281  { vtkMPIController::UseSsendForRMI = (use_send != 0)? 1: 0; }
283 //BTX
284 protected:
286  ~vtkMPIController();
288 
289  // Set the communicator to comm and call InitializeNumberOfProcesses()
290  void InitializeCommunicator(vtkMPICommunicator* comm);
291 
292  // Duplicate the current communicator, creating RMICommunicator
293  void InitializeRMICommunicator();
294 
296 
299  virtual void TriggerRMIInternal(int remoteProcessId,
300  void* arg, int argLength, int rmiTag, bool propagate);
302 
303  // MPI communicator created when Initialize() called.
304  // This is a copy of MPI_COMM_WORLD but uses a new
305  // context, i.e. even if the tags are the same, the
306  // RMI messages will not interfere with user level
307  // messages.
309 
310  friend class vtkMPIOutputWindow;
311 
312  // Initialize only once.
313  static int Initialized;
314 
315  static char ProcessorName[];
316 
318 
319  static int UseSsendForRMI;
320 private:
321  vtkMPIController(const vtkMPIController&); // Not implemented.
322  void operator=(const vtkMPIController&); // Not implemented.
323 //ETX
324 };
326 
327 
328 #endif
329 
330 
int NoBlockSend(const double *data, int length, int remoteProcessId, int tag, vtkMPICommunicator::Request &req)
GLsizei GLsizei GLenum GLenum const GLvoid * data
Definition: vtkgl.h:11339
GLboolean GLuint group
Definition: vtkgl.h:18647
GLsizeiptr size
Definition: vtkgl.h:11843
int NoBlockSend(const unsigned char *data, int length, int remoteProcessId, int tag, vtkMPICommunicator::Request &req)
GLuint GLuint GLsizei GLenum type
Definition: vtkgl.h:11315
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)
GLuint GLuint GLsizei count
Definition: vtkgl.h:11315
GLuint GLsizei GLsizei * length
Definition: vtkgl.h:11992
virtual void TriggerRMIInternal(int remoteProcessId, void *arg, int argLength, int rmiTag, bool propagate)
int vtkIdType
Definition: vtkType.h:268
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:45
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