~ubuntu-branches/ubuntu/trusty/maloc/trusty-proposed

« back to all changes in this revision

Viewing changes to src/psh/vmpi.c

  • Committer: Bazaar Package Importer
  • Author(s): Michael Banck
  • Date: 2006-06-29 15:21:06 UTC
  • Revision ID: james.westby@ubuntu.com-20060629152106-kyqdw6qlc3vmqum3
Tags: upstream-0.2
ImportĀ upstreamĀ versionĀ 0.2

Show diffs side-by-side

added added

removed removed

Lines of Context:
 
1
/*
 
2
 * ***************************************************************************
 
3
 * MALOC = < Minimal Abstraction Layer for Object-oriented C >
 
4
 * Copyright (C) 1994--2006  Michael Holst
 
5
 * 
 
6
 * This program is free software; you can redistribute it and/or modify it
 
7
 * under the terms of the GNU General Public License as published by the
 
8
 * Free Software Foundation; either version 2 of the License, or (at your
 
9
 * option) any later version.
 
10
 * 
 
11
 * This program is distributed in the hope that it will be useful,
 
12
 * but WITHOUT ANY WARRANTY; without even the implied warranty of
 
13
 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
 
14
 * See the GNU General Public License for more details.
 
15
 * 
 
16
 * You should have received a copy of the GNU General Public License along
 
17
 * with this program; if not, write to the Free Software Foundation, Inc.,
 
18
 * 675 Mass Ave, Cambridge, MA 02139, USA.
 
19
 * 
 
20
 * rcsid="$Id: vmpi.c,v 1.20 2006/06/03 07:22:29 mholst Exp $"
 
21
 * ***************************************************************************
 
22
 */
 
23
 
 
24
/*
 
25
 * ***************************************************************************
 
26
 * File:     vmpi.c
 
27
 *
 
28
 * Purpose:  Class Vmpi: methods.
 
29
 *
 
30
 * Notes:    Class Vmpi is a thin object-oriented Clean C layer on top of the
 
31
 *           MPI communication library.  Vmpi provides access to the minimal
 
32
 *           set of ten MPI primitives required to implement the Bank-Holst
 
33
 *           parallel adaptive algorithm, using either the Bank-Holst Oracle
 
34
 *           library, or directly.
 
35
 *
 
36
 *           Vmpi is also just about the smallest possible communication
 
37
 *           library that still contains just enough functionality to solve
 
38
 *           an elliptic equation in parallel using a typical parallel
 
39
 *           algorithm such as domain decomposition or parallel multigrid.
 
40
 *
 
41
 *           This minimal functionality is provided by Vmpi in ten core
 
42
 *           methods: startup, shutdown, size, rank, send, receive,
 
43
 *           broadcast, reduce, barrier, and non-blocking send (sometimes
 
44
 *           very useful).  That is it, and it all fits in about 90 lines
 
45
 *           of object-oriented Clean C code (not counting comments...).
 
46
 *           To fit into such a small amount of code, Vmpi assumes that the
 
47
 *           user does all of his own datapacking into pure character arrays.
 
48
 *           The best way to make use of Vmpi is to pack every message into
 
49
 *           a character buffer in XDR format (RPC's answer to
 
50
 *           machine-independent binary format).  You can accomplish this
 
51
 *           e.g. using the Vio library that sits at the bottom of MC.
 
52
 *
 
53
 *           The 6 primary MPI routines (according to Ian Foster):
 
54
 *           -----------------------------------------------------
 
55
 *
 
56
 *           MPI_Init       : Initiate an MPI computation
 
57
 *           MPI_Finalize   : Terminate a computation
 
58
 *           MPI_Comm_size  : Determine number of processes
 
59
 *           MPI_Comm_rank  : Determine my process identifier
 
60
 *           MPI_Send       : Send a message (blocking)
 
61
 *           MPI_Recv       : Receive a message (blocking)
 
62
 *
 
63
 *           The 4 additional useful MPI routines (according to me):
 
64
 *           -------------------------------------------------------
 
65
 *
 
66
 *           MPI_Barrier    : Synchronizes all processes
 
67
 *           MPI_Bcast      : Sends data from one process to all processes
 
68
 *           MPI_Reduce     : Sums, etc., distributed data (result in root)
 
69
 *           MPI_Isend      : Send a message (non-blocking)
 
70
 *
 
71
 * Author:   Michael Holst
 
72
 * ***************************************************************************
 
73
 */
 
74
 
 
75
#include "vmpi_p.h"
 
76
 
 
77
#if defined(HAVE_MPI_H)
 
78
#   include <mpi.h>
 
79
#endif
 
80
 
 
81
VEMBED(rcsid="$Id: vmpi.c,v 1.20 2006/06/03 07:22:29 mholst Exp $")
 
82
 
 
83
/*
 
84
 * ***************************************************************************
 
85
 * Class Vmpi: Inlineable methods
 
86
 * ***************************************************************************
 
87
 */
 
88
#if !defined(VINLINE_MALOC)
 
89
 
 
90
#endif /* if !defined(VINLINE_MALOC) */
 
91
/*
 
92
 * ***************************************************************************
 
93
 * Class Vmpi: Non-inlineable methods
 
94
 * ***************************************************************************
 
95
 */
 
96
 
 
97
/*
 
98
 * ***************************************************************************
 
99
 * Routine:  Vmpi_init
 
100
 *
 
101
 * Purpose:  The Vmp initializer.
 
102
 *
 
103
 * Author:   Michael Holst
 
104
 * ***************************************************************************
 
105
 */
 
106
VPUBLIC int Vmpi_init(int *argc, char ***argv)
 
107
{
 
108
#if defined(HAVE_MPI_H)
 
109
    return (MPI_SUCCESS == MPI_Init(argc,argv));
 
110
#else
 
111
    return 1;
 
112
#endif
 
113
}
 
114
 
 
115
/*
 
116
 * ***************************************************************************
 
117
 * Routine:  Vmpi_finalize
 
118
 *
 
119
 * Purpose:  The Vmp finalizerr.
 
120
 *
 
121
 * Author:   Michael Holst
 
122
 * ***************************************************************************
 
123
 */
 
124
VPUBLIC int Vmpi_finalize(void)
 
125
{
 
126
#if defined(HAVE_MPI_H)
 
127
    return (MPI_SUCCESS == MPI_Finalize());
 
128
#else
 
129
    return 1;
 
130
#endif
 
131
}
 
132
 
 
133
/*
 
134
 * ***************************************************************************
 
135
 * Routine:  Vmpi_ctor
 
136
 *
 
137
 * Purpose:  The Vmpi constructor.
 
138
 *
 
139
 * Author:   Michael Holst
 
140
 * ***************************************************************************
 
141
 */
 
142
VPUBLIC Vmpi* Vmpi_ctor(void)
 
143
{
 
144
#if defined(HAVE_MPI_H)
 
145
    int dummy;
 
146
#endif
 
147
    Vmpi *thee = VNULL;
 
148
    VDEBUGIO("Vmpi_ctor: CREATING object..");
 
149
    thee = Vmem_malloc( VNULL, 1, sizeof(Vmpi) );
 
150
    thee->mpi_rank = 0;
 
151
    thee->mpi_size = 0;
 
152
#if defined(HAVE_MPI_H)
 
153
    VASSERT( MPI_SUCCESS == MPI_Initialized(&dummy) );
 
154
    VASSERT( MPI_SUCCESS == MPI_Comm_rank(MPI_COMM_WORLD, &(thee->mpi_rank)) );
 
155
    VASSERT( MPI_SUCCESS == MPI_Comm_size(MPI_COMM_WORLD, &(thee->mpi_size)) );
 
156
    Vnm_setIoTag(thee->mpi_rank, thee->mpi_size);
 
157
    Vnm_print(2,"Vmpi_ctor: process %d of %d is ALIVE!\n",
 
158
        thee->mpi_rank, thee->mpi_size);
 
159
#endif
 
160
    VDEBUGIO("..done.\n");
 
161
    return thee;
 
162
}
 
163
 
 
164
/*
 
165
 * ***************************************************************************
 
166
 * Routine:  Vmpi_dtor
 
167
 *
 
168
 * Purpose:  The Vmpi destructor.
 
169
 *
 
170
 * Author:   Michael Holst
 
171
 * ***************************************************************************
 
172
 */
 
173
VPUBLIC void Vmpi_dtor(Vmpi **thee)
 
174
{
 
175
    VDEBUGIO("Vmpi_dtor: DESTROYING object..");
 
176
    Vmem_free( VNULL, 1, sizeof(Vmpi), (void**)thee );
 
177
#if defined(HAVE_MPI_H)
 
178
#if 0
 
179
    VASSERT( MPI_SUCCESS == MPI_Finalize() );
 
180
#endif
 
181
#endif
 
182
    VDEBUGIO("..done.\n");
 
183
}
 
184
 
 
185
/*
 
186
 * ***************************************************************************
 
187
 * Routine:  Vmpi_rank
 
188
 *
 
189
 * Purpose:  Return my processor ID.
 
190
 *
 
191
 * Author:   Michael Holst
 
192
 * ***************************************************************************
 
193
 */
 
194
VPUBLIC int Vmpi_rank(Vmpi *thee)
 
195
{
 
196
    return thee->mpi_rank;
 
197
}
 
198
 
 
199
/*
 
200
 * ***************************************************************************
 
201
 * Routine:  Vmpi_size
 
202
 *
 
203
 * Purpose:  Return the number of processors involved.
 
204
 *
 
205
 * Author:   Michael Holst
 
206
 * ***************************************************************************
 
207
 */
 
208
VPUBLIC int Vmpi_size(Vmpi *thee)
 
209
{
 
210
    return thee->mpi_size;
 
211
}
 
212
 
 
213
/*
 
214
 * ***************************************************************************
 
215
 * Routine:  Vmpi_barr
 
216
 *
 
217
 * Purpose:  An MPI barrier.
 
218
 *
 
219
 * Author:   Michael Holst
 
220
 * ***************************************************************************
 
221
 */
 
222
VPUBLIC int Vmpi_barr(Vmpi *thee)
 
223
{
 
224
    int rc=0;
 
225
#if defined(HAVE_MPI_H)
 
226
    rc = (MPI_SUCCESS == MPI_Barrier(MPI_COMM_WORLD));
 
227
#endif
 
228
    return rc;
 
229
}
 
230
 
 
231
/*
 
232
 * ***************************************************************************
 
233
 * Routine:  Vmpi_send
 
234
 *
 
235
 * Purpose:  An MPI blocking send.
 
236
 *
 
237
 * Author:   Michael Holst
 
238
 * ***************************************************************************
 
239
 */
 
240
VPUBLIC int Vmpi_send(Vmpi *thee, int des, char *buf, int bufsize)
 
241
{
 
242
    int rc=0;
 
243
#if defined(HAVE_MPI_H)
 
244
    int tag=0;
 
245
    rc = (MPI_SUCCESS == MPI_Send(buf,bufsize,MPI_CHAR,des,tag,MPI_COMM_WORLD));
 
246
#endif
 
247
    return rc;
 
248
}
 
249
 
 
250
/*
 
251
 * ***************************************************************************
 
252
 * Routine:  Vmpi_recv
 
253
 *
 
254
 * Purpose:  An MPI blocking receive.
 
255
 *
 
256
 * Author:   Michael Holst
 
257
 * ***************************************************************************
 
258
 */
 
259
VPUBLIC int Vmpi_recv(Vmpi *thee, int src, char *buf, int bufsize)
 
260
{
 
261
    int rc=0;
 
262
#if defined(HAVE_MPI_H)
 
263
    int rsrc=0;
 
264
    int tag=0;
 
265
    MPI_Status stat;
 
266
    if (src < 0) {
 
267
        rsrc = MPI_ANY_SOURCE;
 
268
    } else {
 
269
        rsrc = src;
 
270
    }
 
271
    rc = (MPI_SUCCESS == MPI_Recv(buf,bufsize,MPI_CHAR,rsrc,tag,
 
272
        MPI_COMM_WORLD,&stat));
 
273
#endif
 
274
    return rc;
 
275
}
 
276
 
 
277
/*
 
278
 * ***************************************************************************
 
279
 * Routine:  Vmpi_bcast
 
280
 *
 
281
 * Purpose:  An MPI broadcast.
 
282
 *
 
283
 * Author:   Michael Holst
 
284
 * ***************************************************************************
 
285
 */
 
286
VPUBLIC int Vmpi_bcast(Vmpi *thee, char *buf, int bufsize)
 
287
{
 
288
    int rc=0;
 
289
#if defined(HAVE_MPI_H)
 
290
    rc = (MPI_SUCCESS == MPI_Bcast(buf,bufsize,MPI_CHAR,0,MPI_COMM_WORLD));
 
291
#endif
 
292
    return rc;
 
293
}
 
294
 
 
295
/*
 
296
 * ***************************************************************************
 
297
 * Routine:  Vmpi_reduce
 
298
 *
 
299
 * Purpose:  An MPI reduce.
 
300
 *
 
301
 * Author:   Michael Holst
 
302
 * ***************************************************************************
 
303
 */
 
304
VPUBLIC int Vmpi_reduce(Vmpi *thee, char *sbuf, char *rbuf, int bufsize)
 
305
{
 
306
    int rc=0;
 
307
#if defined(HAVE_MPI_H)
 
308
    rc = (MPI_SUCCESS ==
 
309
        MPI_Reduce(sbuf,rbuf,bufsize,MPI_CHAR,MPI_SUM,0,MPI_COMM_WORLD));
 
310
#endif
 
311
    return rc;
 
312
}
 
313
 
 
314
/*
 
315
 * ***************************************************************************
 
316
 * Routine:  Vmpi_isend
 
317
 *
 
318
 * Purpose:  An MPI non-blocking send.
 
319
 *
 
320
 * Author:   Michael Holst
 
321
 * ***************************************************************************
 
322
 */
 
323
VPUBLIC int Vmpi_isend(Vmpi *thee, int des, char *buf, int bufsize)
 
324
{
 
325
    int rc=0;
 
326
#if defined(HAVE_MPI_H)
 
327
    int tag=0;
 
328
    MPI_Request requ;
 
329
    rc = (MPI_SUCCESS == MPI_Isend(buf,bufsize,MPI_CHAR,des,tag,MPI_COMM_WORLD,
 
330
        &requ));
 
331
#endif
 
332
    return rc;
 
333
}
 
334