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=========================================================================*/
49#ifndef vtkMPIController_h
50#define vtkMPIController_h
51
52#include "vtkParallelMPIModule.h" // For export macro
54// Do not remove this header file. This class contains methods
55// which take arguments defined in vtkMPICommunicator.h by
56// reference.
57#include "vtkMPICommunicator.h" // Needed for direct access to communicator
58
59class vtkIntArray;
60
61class VTKPARALLELMPI_EXPORT vtkMPIController : public vtkMultiProcessController
62{
63
64public:
65
68 void PrintSelf(ostream& os, vtkIndent indent);
69
81 virtual void Initialize(int* argc, char*** argv)
82 { this->Initialize(argc, argv, 0); }
83
84 virtual void Initialize(int* vtkNotUsed(argc), char*** vtkNotUsed(argv),
85 int initializedExternally);
86
90 virtual void Initialize();
91
97 virtual void Finalize() { this->Finalize(0); }
98
99 virtual void Finalize(int finalizedExternally);
100
105 virtual void SingleMethodExecute();
106
112 virtual void MultipleMethodExecute();
113
119 virtual void CreateOutputWindow();
120
125 static char* ErrorString(int err);
126
127
138
140
141 virtual vtkMPIController *PartitionController(int localColor, int localKey);
142
152 int NoBlockSend(const int* data, int length, int remoteProcessId, int tag,
154 { return ((vtkMPICommunicator*)this->Communicator)->NoBlockSend
155 (data ,length, remoteProcessId, tag, req); }
156 int NoBlockSend(const unsigned long* data, int length, int remoteProcessId,
157 int tag, vtkMPICommunicator::Request& req)
158 { return ((vtkMPICommunicator*)this->Communicator)->NoBlockSend
159 (data, length, remoteProcessId, tag, req); }
160 int NoBlockSend(const char* data, int length, int remoteProcessId,
161 int tag, vtkMPICommunicator::Request& req)
162 { return ((vtkMPICommunicator*)this->Communicator)->NoBlockSend
163 (data, length, remoteProcessId, tag, req); }
164 int NoBlockSend( const unsigned char* data, int length, int remoteProcessId,
165 int tag, vtkMPICommunicator::Request& req )
166 { return ((vtkMPICommunicator*)this->Communicator)->NoBlockSend
167 (data, length, remoteProcessId, tag, req);}
168 int NoBlockSend(const float* data, int length, int remoteProcessId,
169 int tag, vtkMPICommunicator::Request& req)
170 { return ((vtkMPICommunicator*)this->Communicator)->NoBlockSend
171 (data, length, remoteProcessId, tag, req); }
172 int NoBlockSend(const double* data, int length, int remoteProcessId,
173 int tag, vtkMPICommunicator::Request& req)
174 { return ((vtkMPICommunicator*)this->Communicator)->NoBlockSend
175 (data, length, remoteProcessId, tag, req); }
176#ifdef VTK_USE_64BIT_IDS
177 int NoBlockSend(const vtkIdType* data, int length, int remoteProcessId,
178 int tag, vtkMPICommunicator::Request& req)
179 { return ((vtkMPICommunicator*)this->Communicator)->NoBlockSend
180 (data, length, remoteProcessId, tag, req); }
181#endif
182
191 int NoBlockReceive(int* data, int length, int remoteProcessId,
192 int tag, vtkMPICommunicator::Request& req)
193 { return ((vtkMPICommunicator*)this->Communicator)->NoBlockReceive
194 (data, length, remoteProcessId, tag, req); }
195 int NoBlockReceive(unsigned long* data, int length,
196 int remoteProcessId, int tag,
198 { return ((vtkMPICommunicator*)this->Communicator)->NoBlockReceive
199 (data, length, remoteProcessId, tag, req); }
200 int NoBlockReceive(char* data, int length, int remoteProcessId,
201 int tag, vtkMPICommunicator::Request& req)
202 { return ((vtkMPICommunicator*)this->Communicator)->NoBlockReceive
203 (data, length, remoteProcessId, tag, req); }
204 int NoBlockReceive(unsigned char* data, int length, int remoteProcessId,
205 int tag, vtkMPICommunicator::Request& req)
206 { return ((vtkMPICommunicator*)this->Communicator)->NoBlockReceive
207 (data, length, remoteProcessId, tag, req); }
208 int NoBlockReceive(float* data, int length, int remoteProcessId,
209 int tag, vtkMPICommunicator::Request& req)
210 { return ((vtkMPICommunicator*)this->Communicator)->NoBlockReceive
211 (data, length, remoteProcessId, tag, req); }
212 int NoBlockReceive(double* data, int length, int remoteProcessId,
213 int tag, vtkMPICommunicator::Request& req)
214 { return ((vtkMPICommunicator*)this->Communicator)->NoBlockReceive
215 (data, length, remoteProcessId, tag, req); }
216#ifdef VTK_USE_64BIT_IDS
217 int NoBlockReceive(vtkIdType* data, int length, int remoteProcessId,
218 int tag, vtkMPICommunicator::Request& req)
219 { return ((vtkMPICommunicator*)this->Communicator)->NoBlockReceive
220 (data, length, remoteProcessId, tag, req); }
221#endif
222
233 int Iprobe(int source, int tag, int* flag, int* actualSource)
234 { return ((vtkMPICommunicator*)this->Communicator)->Iprobe(
235 source, tag, flag, actualSource); }
236 int Iprobe(int source, int tag, int* flag, int* actualSource,
237 int* type, int* size)
238 { return ((vtkMPICommunicator*)this->Communicator)->Iprobe(
239 source, tag, flag, actualSource, type, size); }
240 int Iprobe(int source, int tag, int* flag, int* actualSource,
241 unsigned long* type, int* size)
242 { return ((vtkMPICommunicator*)this->Communicator)->Iprobe(
243 source, tag, flag, actualSource, type, size); }
244 int Iprobe(int source, int tag, int* flag, int* actualSource,
245 const char* type, int* size)
246 { return ((vtkMPICommunicator*)this->Communicator)->Iprobe(
247 source, tag, flag, actualSource, type, size); }
248 int Iprobe(int source, int tag, int* flag, int* actualSource,
249 float* type, int* size)
250 { return ((vtkMPICommunicator*)this->Communicator)->Iprobe(
251 source, tag, flag, actualSource, type, size); }
252 int Iprobe(int source, int tag, int* flag, int* actualSource,
253 double* type, int* size)
254 { return ((vtkMPICommunicator*)this->Communicator)->Iprobe(
255 source, tag, flag, actualSource, type, size); }
256
262 int WaitAll(const int count, vtkMPICommunicator::Request requests[])
263 { return ((vtkMPICommunicator*)this->Communicator)->WaitAll(count,requests);}
264
271 int WaitAny(const int count,vtkMPICommunicator::Request requests[], int& idx)
272 {return ((vtkMPICommunicator*)this->Communicator)->WaitAny(count,requests,idx);}
273
280 const int count, vtkMPICommunicator::Request requests[],
281 vtkIntArray *completed );
282
286 bool TestAll(const int count, vtkMPICommunicator::Request requests[] );
287
294 bool TestAny(const int count,vtkMPICommunicator::Request requests[],int &idx);
295
301 bool TestSome(const int count,vtkMPICommunicator::Request requests[],
302 vtkIntArray *completed );
303
304 static const char* GetProcessorName();
305
310 static void SetUseSsendForRMI(int use_send)
311 { vtkMPIController::UseSsendForRMI = (use_send != 0)? 1: 0; }
313
314protected:
317
318 // Set the communicator to comm and call InitializeNumberOfProcesses()
320
321 // Duplicate the current communicator, creating RMICommunicator
323
329 virtual void TriggerRMIInternal(int remoteProcessId,
330 void* arg, int argLength, int rmiTag, bool propagate);
331
332 // MPI communicator created when Initialize() called.
333 // This is a copy of MPI_COMM_WORLD but uses a new
334 // context, i.e. even if the tags are the same, the
335 // RMI messages will not interfere with user level
336 // messages.
338
339 friend class vtkMPIOutputWindow;
340
341 // Initialize only once.
342 static int Initialized;
343
344 static char ProcessorName[];
345
349 static int UseSsendForRMI;
350
351private:
352 vtkMPIController(const vtkMPIController&) VTK_DELETE_FUNCTION;
353 void operator=(const vtkMPIController&) VTK_DELETE_FUNCTION;
354
355};
356
357
358#endif
359
360
a simple class to control print indentation
Definition: vtkIndent.h:40
dynamic, self-adjusting array of int
Definition: vtkIntArray.h:46
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(int finalizedExternally)
This method is for cleaning up.
virtual void SingleMethodExecute()
Execute the SingleMethod (as define by SetSingleMethod) using this->NumberOfProcesses 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 Initialize(int *vtkNotUsed(argc), char ***vtkNotUsed(argv), int initializedExternally)
This method is for setting up the processes.
void InitializeCommunicator(vtkMPICommunicator *comm)
virtual void Initialize(int *argc, char ***argv)
This method is for setting up the processes.
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)
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 MultipleMethodExecute()
Execute the MultipleMethods (as define by calling SetMultipleMethod for each of the required this->Nu...
static int UseSsendForRMI
When set, TriggerRMI uses Ssend instead of Send.
virtual vtkMPIController * CreateSubController(vtkProcessGroup *group)
Creates a new controller with the processes specified by the given group.
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)
int Iprobe(int source, int tag, int *flag, int *actualSource, float *type, int *size)
virtual vtkMPIController * PartitionController(int localColor, int localKey)
Partitions this controller based on a coloring.
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)
static char * ErrorString(int err)
Given an MPI error code, return a string which contains an error message.
void PrintSelf(ostream &os, vtkIndent indent)
Methods invoked by print to print information about the object including superclasses.
static vtkMPIController * New()
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)
int Iprobe(int source, int tag, int *flag, int *actualSource)
Nonblocking test for a message.
virtual void Finalize()
This method is for cleaning up and has to be called before the end of the program if MPI was initiali...
virtual void TriggerRMIInternal(int remoteProcessId, void *arg, int argLength, int rmiTag, bool propagate)
Implementation for TriggerRMI() provides subclasses an opportunity to modify the behaviour eg.
int NoBlockSend(const double *data, int length, int remoteProcessId, int tag, vtkMPICommunicator::Request &req)
virtual void CreateOutputWindow()
This method can be used to tell the controller to create a special output window in which all message...
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()
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.
virtual void Initialize(int *vtkNotUsed(argc), char ***vtkNotUsed(argv))=0
This method is for setting up the processes.
A subgroup of processes from a communicator.
@ length
Definition: vtkX3D.h:393
@ type
Definition: vtkX3D.h:516
@ size
Definition: vtkX3D.h:253
@ data
Definition: vtkX3D.h:315
boost::graph_traits< vtkGraph * >::vertex_descriptor source(boost::graph_traits< vtkGraph * >::edge_descriptor e, vtkGraph *)
int vtkIdType
Definition: vtkType.h:287