Statistiques
| Révision :

root / src / grid / HPL_reduce.c @ 1

Historique | Voir | Annoter | Télécharger (7,31 ko)

1
/* 
2
 * -- High Performance Computing Linpack Benchmark (HPL)                
3
 *    HPL - 2.0 - September 10, 2008                          
4
 *    Antoine P. Petitet                                                
5
 *    University of Tennessee, Knoxville                                
6
 *    Innovative Computing Laboratory                                 
7
 *    (C) Copyright 2000-2008 All Rights Reserved                       
8
 *                                                                      
9
 * -- Copyright notice and Licensing terms:                             
10
 *                                                                      
11
 * Redistribution  and  use in  source and binary forms, with or without
12
 * modification, are  permitted provided  that the following  conditions
13
 * are met:                                                             
14
 *                                                                      
15
 * 1. Redistributions  of  source  code  must retain the above copyright
16
 * notice, this list of conditions and the following disclaimer.        
17
 *                                                                      
18
 * 2. Redistributions in binary form must reproduce  the above copyright
19
 * notice, this list of conditions,  and the following disclaimer in the
20
 * documentation and/or other materials provided with the distribution. 
21
 *                                                                      
22
 * 3. All  advertising  materials  mentioning  features  or  use of this
23
 * software must display the following acknowledgement:                 
24
 * This  product  includes  software  developed  at  the  University  of
25
 * Tennessee, Knoxville, Innovative Computing Laboratory.             
26
 *                                                                      
27
 * 4. The name of the  University,  the name of the  Laboratory,  or the
28
 * names  of  its  contributors  may  not  be used to endorse or promote
29
 * products  derived   from   this  software  without  specific  written
30
 * permission.                                                          
31
 *                                                                      
32
 * -- Disclaimer:                                                       
33
 *                                                                      
34
 * THIS  SOFTWARE  IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
35
 * ``AS IS'' AND ANY EXPRESS OR IMPLIED WARRANTIES,  INCLUDING,  BUT NOT
36
 * LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR
37
 * A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE UNIVERSITY
38
 * OR  CONTRIBUTORS  BE  LIABLE FOR ANY  DIRECT,  INDIRECT,  INCIDENTAL,
39
 * SPECIAL,  EXEMPLARY,  OR  CONSEQUENTIAL DAMAGES  (INCLUDING,  BUT NOT
40
 * LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE,
41
 * DATA OR PROFITS; OR BUSINESS INTERRUPTION)  HOWEVER CAUSED AND ON ANY
42
 * THEORY OF LIABILITY, WHETHER IN CONTRACT,  STRICT LIABILITY,  OR TORT
43
 * (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE
44
 * OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE. 
45
 * ---------------------------------------------------------------------
46
 */ 
47
/*
48
 * Include files
49
 */
50
#include "hpl.h"
51

    
52
#ifdef STDC_HEADERS
53
int HPL_reduce
54
(
55
   void *                           BUFFER,
56
   const int                        COUNT,
57
   const HPL_T_TYPE                 DTYPE,
58
   const HPL_T_OP                   OP,
59
   const int                        ROOT,
60
   MPI_Comm                         COMM
61
)
62
#else
63
int HPL_reduce
64
( BUFFER, COUNT, DTYPE, OP, ROOT, COMM )
65
   void *                           BUFFER;
66
   const int                        COUNT;
67
   const HPL_T_TYPE                 DTYPE;
68
   const HPL_T_OP                   OP;
69
   const int                        ROOT;
70
   MPI_Comm                         COMM;
71
#endif
72
{
73
/* 
74
 * Purpose
75
 * =======
76
 *
77
 * HPL_reduce performs a global reduce operation across all processes of
78
 * a group.  Note that the input buffer is  used as workarray and in all
79
 * processes but the accumulating process corrupting the original data.
80
 *
81
 * Arguments
82
 * =========
83
 *
84
 * BUFFER  (local input/output)          void *
85
 *         On entry,  BUFFER  points to  the  buffer to be  reduced.  On
86
 *         exit,  and  in process of rank  ROOT  this array contains the
87
 *         reduced data.  This  buffer  is also used as workspace during
88
 *         the operation in the other processes of the group.
89
 *
90
 * COUNT   (global input)                const int
91
 *         On entry,  COUNT  indicates the number of entries in  BUFFER.
92
 *         COUNT must be at least zero.
93
 *
94
 * DTYPE   (global input)                const HPL_T_TYPE
95
 *         On entry,  DTYPE  specifies the type of the buffers operands.
96
 *
97
 * OP      (global input)                const HPL_T_OP 
98
 *         On entry, OP is a pointer to the local combine function.
99
 *
100
 * ROOT    (global input)                const int
101
 *         On entry, ROOT is the coordinate of the accumulating process.
102
 *
103
 * COMM    (global/local input)          MPI_Comm
104
 *         The MPI communicator identifying the process collection.
105
 *
106
 * ---------------------------------------------------------------------
107
 */ 
108
/*
109
 * .. Local Variables ..
110
 */
111
   MPI_Status                 status;
112
   void                       * buffer = NULL;
113
   int                        hplerr=MPI_SUCCESS, d=1, i, ip2=1, mask=0,
114
                              mpierr, mydist, partner, rank, size, 
115
                              tag = MSGID_BEGIN_COLL;
116
/* ..
117
 * .. Executable Statements ..
118
 */
119
   if( COUNT <= 0 ) return( MPI_SUCCESS );
120
   mpierr = MPI_Comm_size( COMM, &size );
121
   if( size  == 1 ) return( MPI_SUCCESS );
122
   mpierr = MPI_Comm_rank( COMM, &rank );
123
   i = size - 1; while( i > 1 ) { i >>= 1; d++; }
124

    
125
   if( DTYPE == HPL_INT )
126
      buffer = (void *)( (int *)   malloc( (size_t)(COUNT) * 
127
                                           sizeof( int    ) ) );
128
   else
129
      buffer = (void *)( (double *)malloc( (size_t)(COUNT) *
130
                                           sizeof( double ) ) );
131

    
132
   if( !( buffer ) )
133
   { HPL_pabort( __LINE__, "HPL_reduce", "Memory allocation failed" ); }
134

    
135
   if( ( mydist = MModSub( rank, ROOT, size ) ) == 0 )
136
   {
137
      do
138
      {
139
         mpierr = MPI_Recv( buffer, COUNT, HPL_2_MPI_TYPE( DTYPE ),
140
                            MModAdd( ROOT, ip2, size ), tag, COMM,
141
                            &status );
142
         if( mpierr != MPI_SUCCESS ) hplerr = mpierr;
143
         OP( COUNT, buffer, BUFFER, DTYPE );
144
         ip2 <<= 1; d--;
145
      } while( d );
146
   }
147
   else
148
   {
149
      do
150
      {
151
         if( ( mydist & mask ) == 0 )
152
         {
153
            partner = mydist ^ ip2;
154

    
155
            if( mydist & ip2 )
156
            {
157
               partner = MModAdd( ROOT, partner, size );
158
               mpierr = MPI_Send( BUFFER, COUNT, HPL_2_MPI_TYPE( DTYPE ),
159
                                  partner, tag, COMM );
160
            }
161
            else if( partner < size )
162
            {
163
               partner = MModAdd( ROOT, partner, size );
164
               mpierr  = MPI_Recv( buffer, COUNT, HPL_2_MPI_TYPE( DTYPE ),
165
                                   partner, tag, COMM, &status );
166
               OP( COUNT, buffer, BUFFER, DTYPE );
167
            }
168
            if( mpierr != MPI_SUCCESS ) hplerr = mpierr;
169
         }
170
         mask ^= ip2; ip2 <<= 1; d--;
171
      } while( d );
172
   }
173
   if( buffer ) free( buffer );
174

    
175
   return( hplerr );
176
/*
177
 * End of HPL_reduce
178
 */
179
}