Statistiques
| Révision :

root / src / grid / HPL_broadcast.c @ 9

Historique | Voir | Annoter | Télécharger (6 ko)

1 1 equemene
/*
2 1 equemene
 * -- High Performance Computing Linpack Benchmark (HPL)
3 1 equemene
 *    HPL - 2.0 - September 10, 2008
4 1 equemene
 *    Antoine P. Petitet
5 1 equemene
 *    University of Tennessee, Knoxville
6 1 equemene
 *    Innovative Computing Laboratory
7 1 equemene
 *    (C) Copyright 2000-2008 All Rights Reserved
8 1 equemene
 *
9 1 equemene
 * -- Copyright notice and Licensing terms:
10 1 equemene
 *
11 1 equemene
 * Redistribution  and  use in  source and binary forms, with or without
12 1 equemene
 * modification, are  permitted provided  that the following  conditions
13 1 equemene
 * are met:
14 1 equemene
 *
15 1 equemene
 * 1. Redistributions  of  source  code  must retain the above copyright
16 1 equemene
 * notice, this list of conditions and the following disclaimer.
17 1 equemene
 *
18 1 equemene
 * 2. Redistributions in binary form must reproduce  the above copyright
19 1 equemene
 * notice, this list of conditions,  and the following disclaimer in the
20 1 equemene
 * documentation and/or other materials provided with the distribution.
21 1 equemene
 *
22 1 equemene
 * 3. All  advertising  materials  mentioning  features  or  use of this
23 1 equemene
 * software must display the following acknowledgement:
24 1 equemene
 * This  product  includes  software  developed  at  the  University  of
25 1 equemene
 * Tennessee, Knoxville, Innovative Computing Laboratory.
26 1 equemene
 *
27 1 equemene
 * 4. The name of the  University,  the name of the  Laboratory,  or the
28 1 equemene
 * names  of  its  contributors  may  not  be used to endorse or promote
29 1 equemene
 * products  derived   from   this  software  without  specific  written
30 1 equemene
 * permission.
31 1 equemene
 *
32 1 equemene
 * -- Disclaimer:
33 1 equemene
 *
34 1 equemene
 * THIS  SOFTWARE  IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
35 1 equemene
 * ``AS IS'' AND ANY EXPRESS OR IMPLIED WARRANTIES,  INCLUDING,  BUT NOT
36 1 equemene
 * LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR
37 1 equemene
 * A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE UNIVERSITY
38 1 equemene
 * OR  CONTRIBUTORS  BE  LIABLE FOR ANY  DIRECT,  INDIRECT,  INCIDENTAL,
39 1 equemene
 * SPECIAL,  EXEMPLARY,  OR  CONSEQUENTIAL DAMAGES  (INCLUDING,  BUT NOT
40 1 equemene
 * LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE,
41 1 equemene
 * DATA OR PROFITS; OR BUSINESS INTERRUPTION)  HOWEVER CAUSED AND ON ANY
42 1 equemene
 * THEORY OF LIABILITY, WHETHER IN CONTRACT,  STRICT LIABILITY,  OR TORT
43 1 equemene
 * (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE
44 1 equemene
 * OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
45 1 equemene
 * ---------------------------------------------------------------------
46 1 equemene
 */
47 1 equemene
/*
48 1 equemene
 * Include files
49 1 equemene
 */
50 1 equemene
#include "hpl.h"
51 1 equemene
52 1 equemene
#ifdef STDC_HEADERS
53 1 equemene
int HPL_broadcast
54 1 equemene
(
55 1 equemene
   void *                           BUFFER,
56 1 equemene
   const int                        COUNT,
57 1 equemene
   const HPL_T_TYPE                 DTYPE,
58 1 equemene
   const int                        ROOT,
59 1 equemene
   MPI_Comm                         COMM
60 1 equemene
)
61 1 equemene
#else
62 1 equemene
int HPL_broadcast
63 1 equemene
( BUFFER, COUNT, DTYPE, ROOT, COMM )
64 1 equemene
   void *                           BUFFER;
65 1 equemene
   const int                        COUNT;
66 1 equemene
   const HPL_T_TYPE                 DTYPE;
67 1 equemene
   const int                        ROOT;
68 1 equemene
   MPI_Comm                         COMM;
69 1 equemene
#endif
70 1 equemene
{
71 1 equemene
/*
72 1 equemene
 * Purpose
73 1 equemene
 * =======
74 1 equemene
 *
75 1 equemene
 * HPL_broadcast broadcasts  a message from the process with rank ROOT to
76 1 equemene
 * all processes in the group.
77 1 equemene
 *
78 1 equemene
 * Arguments
79 1 equemene
 * =========
80 1 equemene
 *
81 1 equemene
 * BUFFER  (local input/output)          void *
82 1 equemene
 *         On entry,  BUFFER  points to  the  buffer to be broadcast. On
83 1 equemene
 *         exit, this array contains the broadcast data and is identical
84 1 equemene
 *         on all processes in the group.
85 1 equemene
 *
86 1 equemene
 * COUNT   (global input)                const int
87 1 equemene
 *         On entry,  COUNT  indicates the number of entries in  BUFFER.
88 1 equemene
 *         COUNT must be at least zero.
89 1 equemene
 *
90 1 equemene
 * DTYPE   (global input)                const HPL_T_TYPE
91 1 equemene
 *         On entry,  DTYPE  specifies the type of the buffers operands.
92 1 equemene
 *
93 1 equemene
 * ROOT    (global input)                const int
94 1 equemene
 *         On entry, ROOT is the coordinate of the source process.
95 1 equemene
 *
96 1 equemene
 * COMM    (global/local input)          MPI_Comm
97 1 equemene
 *         The MPI communicator identifying the process collection.
98 1 equemene
 *
99 1 equemene
 * ---------------------------------------------------------------------
100 1 equemene
 */
101 1 equemene
/*
102 1 equemene
 * .. Local Variables ..
103 1 equemene
 */
104 1 equemene
   int                        hplerr=MPI_SUCCESS, ip2=1, kk, mask=1,
105 1 equemene
                              mpierr, mydist, partner, rank, size,
106 1 equemene
                              tag = MSGID_BEGIN_COLL;
107 1 equemene
   MPI_Status                 status;
108 1 equemene
/* ..
109 1 equemene
 * .. Executable Statements ..
110 1 equemene
 */
111 1 equemene
   if( COUNT <= 0 ) return( MPI_SUCCESS );
112 1 equemene
   mpierr = MPI_Comm_size( COMM, &size ); if( size <= 1 ) return( mpierr );
113 1 equemene
   mpierr = MPI_Comm_rank( COMM, &rank );
114 1 equemene
115 1 equemene
   kk = size - 1;
116 1 equemene
   while( kk > 1 ) { kk >>= 1; ip2 <<= 1; mask <<= 1; mask++; }
117 1 equemene
   mydist = MModSub( rank, ROOT, size );
118 1 equemene
119 1 equemene
   do
120 1 equemene
   {
121 1 equemene
      mask ^= ip2;
122 1 equemene
      if( ( mydist & mask ) == 0 )
123 1 equemene
      {
124 1 equemene
         partner = mydist ^ ip2;
125 1 equemene
126 1 equemene
         if( mydist & ip2 )
127 1 equemene
         {
128 1 equemene
            partner = MModAdd( ROOT, partner, size );
129 1 equemene
            mpierr  = MPI_Recv(  BUFFER, COUNT, HPL_2_MPI_TYPE( DTYPE ),
130 1 equemene
                                 partner, tag, COMM, &status );
131 1 equemene
         }
132 1 equemene
         else if( partner < size )
133 1 equemene
         {
134 1 equemene
            partner = MModAdd( ROOT, partner, size );
135 1 equemene
            mpierr  = MPI_Send( BUFFER, COUNT, HPL_2_MPI_TYPE( DTYPE ),
136 1 equemene
                                partner, tag, COMM );
137 1 equemene
         }
138 1 equemene
         if( mpierr != MPI_SUCCESS ) hplerr = mpierr;
139 1 equemene
      }
140 1 equemene
      ip2 >>= 1;
141 1 equemene
   } while( ip2 );
142 1 equemene
143 1 equemene
   return( hplerr );
144 1 equemene
/*
145 1 equemene
 * End of HPL_broadcast
146 1 equemene
 */
147 1 equemene
}