VTK  9.0.2
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 =========================================================================*/
43 #ifndef vtkMPIController_h
44 #define vtkMPIController_h
45 
47 #include "vtkParallelMPIModule.h" // For export macro
48 // Do not remove this header file. This class contains methods
49 // which take arguments defined in vtkMPICommunicator.h by
50 // reference.
51 #include "vtkMPICommunicator.h" // Needed for direct access to communicator
52 
53 class vtkIntArray;
54 
55 class VTKPARALLELMPI_EXPORT vtkMPIController : public vtkMultiProcessController
56 {
57 
58 public:
59  static vtkMPIController* New();
61  void PrintSelf(ostream& os, vtkIndent indent) override;
62 
74  virtual void Initialize(int* argc, char*** argv) override { this->Initialize(argc, argv, 0); }
75 
76  virtual void Initialize(
77  int* vtkNotUsed(argc), char*** vtkNotUsed(argv), int initializedExternally) override;
78 
82  virtual void Initialize();
83 
89  virtual void Finalize() override { this->Finalize(0); }
90 
91  virtual void Finalize(int finalizedExternally) override;
92 
97  virtual void SingleMethodExecute() override;
98 
104  virtual void MultipleMethodExecute() override;
105 
111  virtual void CreateOutputWindow() override;
112 
117  static char* ErrorString(int err);
118 
129 
131 
132  virtual vtkMPIController* PartitionController(int localColor, int localKey) override;
133 
144  const int* data, int length, int remoteProcessId, int tag, vtkMPICommunicator::Request& req)
145  {
146  return ((vtkMPICommunicator*)this->Communicator)
147  ->NoBlockSend(data, length, remoteProcessId, tag, req);
148  }
149  int NoBlockSend(const unsigned long* data, int length, int remoteProcessId, int tag,
151  {
152  return ((vtkMPICommunicator*)this->Communicator)
153  ->NoBlockSend(data, length, remoteProcessId, tag, req);
154  }
156  const char* data, int length, int remoteProcessId, int tag, vtkMPICommunicator::Request& req)
157  {
158  return ((vtkMPICommunicator*)this->Communicator)
159  ->NoBlockSend(data, length, remoteProcessId, tag, req);
160  }
161  int NoBlockSend(const unsigned char* data, int length, int remoteProcessId, int tag,
163  {
164  return ((vtkMPICommunicator*)this->Communicator)
165  ->NoBlockSend(data, length, remoteProcessId, tag, req);
166  }
168  const float* data, int length, int remoteProcessId, int tag, vtkMPICommunicator::Request& req)
169  {
170  return ((vtkMPICommunicator*)this->Communicator)
171  ->NoBlockSend(data, length, remoteProcessId, tag, req);
172  }
174  const double* data, int length, int remoteProcessId, int tag, vtkMPICommunicator::Request& req)
175  {
176  return ((vtkMPICommunicator*)this->Communicator)
177  ->NoBlockSend(data, length, remoteProcessId, tag, req);
178  }
179 #ifdef VTK_USE_64BIT_IDS
180  int NoBlockSend(const vtkIdType* data, int length, int remoteProcessId, int tag,
182  {
183  return ((vtkMPICommunicator*)this->Communicator)
184  ->NoBlockSend(data, length, remoteProcessId, tag, req);
185  }
186 #endif
187 
197  int* data, int length, int remoteProcessId, int tag, vtkMPICommunicator::Request& req)
198  {
199  return ((vtkMPICommunicator*)this->Communicator)
200  ->NoBlockReceive(data, length, remoteProcessId, tag, req);
201  }
203  unsigned long* data, int length, int remoteProcessId, int tag, vtkMPICommunicator::Request& req)
204  {
205  return ((vtkMPICommunicator*)this->Communicator)
206  ->NoBlockReceive(data, length, remoteProcessId, tag, req);
207  }
209  char* data, int length, int remoteProcessId, int tag, vtkMPICommunicator::Request& req)
210  {
211  return ((vtkMPICommunicator*)this->Communicator)
212  ->NoBlockReceive(data, length, remoteProcessId, tag, req);
213  }
215  unsigned char* data, int length, int remoteProcessId, int tag, vtkMPICommunicator::Request& req)
216  {
217  return ((vtkMPICommunicator*)this->Communicator)
218  ->NoBlockReceive(data, length, remoteProcessId, tag, req);
219  }
221  float* data, int length, int remoteProcessId, int tag, vtkMPICommunicator::Request& req)
222  {
223  return ((vtkMPICommunicator*)this->Communicator)
224  ->NoBlockReceive(data, length, remoteProcessId, tag, req);
225  }
227  double* data, int length, int remoteProcessId, int tag, vtkMPICommunicator::Request& req)
228  {
229  return ((vtkMPICommunicator*)this->Communicator)
230  ->NoBlockReceive(data, length, remoteProcessId, tag, req);
231  }
232 #ifdef VTK_USE_64BIT_IDS
233  int NoBlockReceive(
234  vtkIdType* data, int length, int remoteProcessId, int tag, vtkMPICommunicator::Request& req)
235  {
236  return ((vtkMPICommunicator*)this->Communicator)
237  ->NoBlockReceive(data, length, remoteProcessId, tag, req);
238  }
239 #endif
240 
251  int Iprobe(int source, int tag, int* flag, int* actualSource)
252  {
253  return ((vtkMPICommunicator*)this->Communicator)->Iprobe(source, tag, flag, actualSource);
254  }
255  int Iprobe(int source, int tag, int* flag, int* actualSource, int* type, int* size)
256  {
257  return ((vtkMPICommunicator*)this->Communicator)
258  ->Iprobe(source, tag, flag, actualSource, type, size);
259  }
260  int Iprobe(int source, int tag, int* flag, int* actualSource, unsigned long* type, int* size)
261  {
262  return ((vtkMPICommunicator*)this->Communicator)
263  ->Iprobe(source, tag, flag, actualSource, type, size);
264  }
265  int Iprobe(int source, int tag, int* flag, int* actualSource, const char* type, int* size)
266  {
267  return ((vtkMPICommunicator*)this->Communicator)
268  ->Iprobe(source, tag, flag, actualSource, type, size);
269  }
270  int Iprobe(int source, int tag, int* flag, int* actualSource, float* type, int* size)
271  {
272  return ((vtkMPICommunicator*)this->Communicator)
273  ->Iprobe(source, tag, flag, actualSource, type, size);
274  }
275  int Iprobe(int source, int tag, int* flag, int* actualSource, double* type, int* size)
276  {
277  return ((vtkMPICommunicator*)this->Communicator)
278  ->Iprobe(source, tag, flag, actualSource, type, size);
279  }
280 
286  int WaitAll(const int count, vtkMPICommunicator::Request requests[])
287  {
288  return ((vtkMPICommunicator*)this->Communicator)->WaitAll(count, requests);
289  }
290 
297  int WaitAny(const int count, vtkMPICommunicator::Request requests[], int& idx)
298  VTK_SIZEHINT(requests, count)
299  {
300  return ((vtkMPICommunicator*)this->Communicator)->WaitAny(count, requests, idx);
301  }
302 
308  int WaitSome(const int count, vtkMPICommunicator::Request requests[], vtkIntArray* completed)
309  VTK_SIZEHINT(requests, count);
310 
314  bool TestAll(const int count, vtkMPICommunicator::Request requests[]);
315 
322  bool TestAny(const int count, vtkMPICommunicator::Request requests[], int& idx)
323  VTK_SIZEHINT(requests, count);
324 
330  bool TestSome(const int count, vtkMPICommunicator::Request requests[], vtkIntArray* completed)
331  VTK_SIZEHINT(requests, count);
332 
333  static const char* GetProcessorName();
334 
339  static void SetUseSsendForRMI(int use_send)
340  {
341  vtkMPIController::UseSsendForRMI = (use_send != 0) ? 1 : 0;
342  }
344 
345 protected:
347  ~vtkMPIController() override;
348 
349  // Set the communicator to comm and call InitializeNumberOfProcesses()
351 
352  // Duplicate the current communicator, creating RMICommunicator
354 
360  virtual void TriggerRMIInternal(
361  int remoteProcessId, void* arg, int argLength, int rmiTag, bool propagate) override;
362 
363  // MPI communicator created when Initialize() called.
364  // This is a copy of MPI_COMM_WORLD but uses a new
365  // context, i.e. even if the tags are the same, the
366  // RMI messages will not interfere with user level
367  // messages.
369 
370  friend class vtkMPIOutputWindow;
371 
372  // Initialize only once.
373  static int Initialized;
374 
375  static char ProcessorName[];
376 
380  static int UseSsendForRMI;
381 
382 private:
383  vtkMPIController(const vtkMPIController&) = delete;
384  void operator=(const vtkMPIController&) = delete;
385 };
386 
387 #endif
a simple class to control print indentation
Definition: vtkIndent.h:34
dynamic, self-adjusting array of int
Definition: vtkIntArray.h:40
Class for creating user defined MPI communicators.
Process communication using MPI.
int Iprobe(int source, int tag, int *flag, int *actualSource, int *type, int *size)
int WaitAll(const int count, vtkMPICommunicator::Request requests[])
Given the request objects of a set of non-blocking operations (send and/or receive) this method block...
int NoBlockReceive(char *data, int length, int remoteProcessId, int tag, vtkMPICommunicator::Request &req)
virtual void Finalize() override
This method is for cleaning up and has to be called before the end of the program if MPI was initiali...
~vtkMPIController() override
virtual void Initialize(int *vtkNotUsed(argc), char ***vtkNotUsed(argv), int initializedExternally) override
This method is for setting up the processes.
int WaitSome(const int count, vtkMPICommunicator::Request requests[], vtkIntArray *completed)
Blocks until one or more of the specified requests in the given request request array completes.
virtual void Finalize(int finalizedExternally) override
This method is for cleaning up.
virtual void TriggerRMIInternal(int remoteProcessId, void *arg, int argLength, int rmiTag, bool propagate) override
Implementation for TriggerRMI() provides subclasses an opportunity to modify the behaviour eg.
void InitializeCommunicator(vtkMPICommunicator *comm)
int NoBlockSend(const char *data, int length, int remoteProcessId, int tag, vtkMPICommunicator::Request &req)
int Iprobe(int source, int tag, int *flag, int *actualSource, double *type, int *size)
int NoBlockSend(const unsigned char *data, int length, int remoteProcessId, int tag, vtkMPICommunicator::Request &req)
virtual vtkMPIController * CreateSubController(vtkProcessGroup *group) override
Creates a new controller with the processes specified by the given group.
int Iprobe(int source, int tag, int *flag, int *actualSource, const char *type, int *size)
int NoBlockReceive(float *data, int length, int remoteProcessId, int tag, vtkMPICommunicator::Request &req)
virtual void CreateOutputWindow() override
This method can be used to tell the controller to create a special output window in which all message...
static vtkMPIController * New()
static int UseSsendForRMI
When set, TriggerRMI uses Ssend instead of Send.
bool TestSome(const int count, vtkMPICommunicator::Request requests[], vtkIntArray *completed)
Return true iff one or more of the communicator request objects is complete.
bool TestAny(const int count, vtkMPICommunicator::Request requests[], int &idx)
Returns true iff at least one of the communication request objects is complete.
static vtkMPICommunicator * WorldRMICommunicator
void SetCommunicator(vtkMPICommunicator *comm)
MPIController uses this communicator in all sends and receives.
int NoBlockSend(const float *data, int length, int remoteProcessId, int tag, vtkMPICommunicator::Request &req)
virtual void Initialize(int *argc, char ***argv) override
This method is for setting up the processes.
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
int Iprobe(int source, int tag, int *flag, int *actualSource, float *type, int *size)
int NoBlockReceive(double *data, int length, int remoteProcessId, int tag, vtkMPICommunicator::Request &req)
static const char * GetProcessorName()
int NoBlockReceive(unsigned long *data, int length, int remoteProcessId, int tag, vtkMPICommunicator::Request &req)
virtual void MultipleMethodExecute() override
Execute the MultipleMethods (as define by calling SetMultipleMethod for each of the required this->Nu...
int NoBlockSend(const int *data, int length, int remoteProcessId, int tag, vtkMPICommunicator::Request &req)
This method sends data to another process (non-blocking).
int WaitAny(const int count, vtkMPICommunicator::Request requests[], int &idx)
Blocks until one of the specified requests in the given request array completes.
static void SetUseSsendForRMI(int use_send)
When set to 1, TriggerRMI uses Ssend() instead of Send() calls.
int NoBlockReceive(int *data, int length, int remoteProcessId, int tag, vtkMPICommunicator::Request &req)
This method receives data from a corresponding send (non-blocking).
int NoBlockReceive(unsigned char *data, int length, int remoteProcessId, int tag, vtkMPICommunicator::Request &req)
virtual vtkMPIController * PartitionController(int localColor, int localKey) override
Partitions this controller based on a coloring.
int Iprobe(int source, int tag, int *flag, int *actualSource)
Nonblocking test for a message.
int NoBlockSend(const double *data, int length, int remoteProcessId, int tag, vtkMPICommunicator::Request &req)
virtual void SingleMethodExecute() override
Execute the SingleMethod (as define by SetSingleMethod) using this->NumberOfProcesses processes.
bool TestAll(const int count, vtkMPICommunicator::Request requests[])
Returns true iff all of the communication request objects are complete.
static int Initialized
virtual void Initialize()
Same as Initialize(0, 0, 1).
static int GetUseSsendForRMI()
static char * ErrorString(int err)
Given an MPI error code, return a string which contains an error message.
int NoBlockSend(const unsigned long *data, int length, int remoteProcessId, int tag, vtkMPICommunicator::Request &req)
int Iprobe(int source, int tag, int *flag, int *actualSource, unsigned long *type, int *size)
void InitializeRMICommunicator()
Multiprocessing communication superclass.
A subgroup of processes from a communicator.
@ length
Definition: vtkX3D.h:399
@ type
Definition: vtkX3D.h:522
@ size
Definition: vtkX3D.h:259
@ data
Definition: vtkX3D.h:321
boost::graph_traits< vtkGraph * >::vertex_descriptor source(boost::graph_traits< vtkGraph * >::edge_descriptor e, vtkGraph *)
int vtkIdType
Definition: vtkType.h:338
#define VTK_SIZEHINT(...)