Statistiques
| Révision :

root / src / pauxil / HPL_pdlaprnt.c @ 1

Historique | Voir | Annoter | Télécharger (7,99 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
void HPL_pdlaprnt
54
(
55
   const HPL_T_grid *               GRID,
56
   const int                        M,
57
   const int                        N,
58
   const int                        NB,
59
   double *                         A,
60
   const int                        LDA,
61
   const int                        IAROW,
62
   const int                        IACOL,
63
   const char *                     CMATNM
64
)
65
#else
66
void HPL_pdlaprnt
67
( GRID, M, N, NB, A, LDA, IAROW, IACOL, CMATNM )
68
   const HPL_T_grid *               GRID;
69
   const int                        M;
70
   const int                        N;
71
   const int                        NB;
72
   double *                         A;
73
   const int                        LDA;
74
   const int                        IAROW;
75
   const int                        IACOL;
76
   const char *                     CMATNM;
77
#endif
78
{
79
/* 
80
 * Purpose
81
 * =======
82
 *
83
 * HPL_pdlaprnt prints  to  standard  error a distributed matrix A. The
84
 * local pieces of  A  are sent to the process of coordinates  (0,0)  in
85
 * the grid and then printed.
86
 *
87
 * Arguments
88
 * =========
89
 *
90
 * GRID    (local input)                 const HPL_T_grid *
91
 *         On entry,  GRID  points  to the data structure containing the
92
 *         process grid information.
93
 *
94
 * M       (global input)                const int
95
 *         On entry,  M  specifies the number of rows of the coefficient
96
 *         matrix A. M must be at least zero.
97
 *
98
 * N       (global input)                const int
99
 *         On  entry,   N   specifies  the  number  of  columns  of  the
100
 *         coefficient matrix A. N must be at least zero.
101
 *
102
 * NB      (global input)                const int
103
 *         On entry,  NB specifies the blocking factor used to partition
104
 *         and distribute the matrix. NB must be larger than one.
105
 *
106
 * A       (local input)                 double *
107
 *         On entry,  A  points to an  array of dimension (LDA,LocQ(N)).
108
 *         This array contains the coefficient matrix to be printed.
109
 *
110
 * LDA     (local input)                 const int
111
 *         On entry, LDA specifies the leading dimension of the array A.
112
 *         LDA must be at least max(1,LocP(M)).
113
 *
114
 * IAROW   (global input)                const int
115
 *         On entry,  IAROW  specifies the row process coordinate owning
116
 *         the  first row of A.  IAROW  must be  larger than or equal to
117
 *         zero and less than NPROW.
118
 *
119
 * IACOL   (global input)                const int
120
 *         On entry,  IACOL  specifies  the  column  process  coordinate
121
 *         owning the  first column  of A. IACOL  must be larger than or
122
 *         equal to zero and less than NPCOL.
123
 *
124
 * CMATNM  (global input)                const char *
125
 *         On entry, CMATNM is the name of the matrix to be printed.
126
 *
127
 * ---------------------------------------------------------------------
128
 */ 
129
/*
130
 * .. Local Variables ..
131
 */
132
   MPI_Comm                   Acomm;
133
   double                     * buf = NULL;
134
   int                        h, i, ib, icurcol=IACOL, icurrow=IAROW,
135
                              ii=0, j, jb, jj=0, mycol, myrow, npcol,
136
                              nprow, src;
137
/* ..
138
 * .. Executable Statements ..
139
 */
140
   (void) HPL_grid_info( GRID, &nprow, &npcol, &myrow, &mycol );
141
   Acomm = GRID->all_comm; 
142
   if( ( myrow == 0 ) && ( mycol == 0 ) )
143
      buf = (double*)malloc( (size_t)(NB) * sizeof( double ) );
144

    
145
   for( j = 0; j < N; j += NB )
146
   {
147
      jb = N-j; jb = Mmin( jb, NB );
148
      for( h = 0; h < jb; h++ )
149
      {
150
         (void) HPL_barrier( Acomm );
151

    
152
         for( i = 0; i < M; i += NB )
153
         {
154
            ib = M-i; ib = Mmin( ib, NB );
155
            if( ( icurrow == 0 ) && ( icurcol == 0 ) )
156
            {
157
               if( ( myrow == 0 ) && ( mycol == 0 ) )
158
                  HPL_dlaprnt( ib, 1, Mptr( A, ii, jj+h, LDA ), i+1,
159
                               j+h+1, LDA, CMATNM );
160
            }
161
            else
162
            {
163
               if( ( myrow == icurrow ) && ( mycol == icurcol ) )
164
               {
165
                  (void) HPL_send( Mptr( A, ii, jj+h, LDA ), ib, 0,
166
                                   9000+(j+h)*M+i, Acomm );
167
               }
168
               else if( ( myrow == 0 ) && ( mycol == 0 ) )
169
               {
170
                  src = HPL_pnum( GRID, icurrow, icurcol );
171
                  (void) HPL_recv( buf, ib, src, 9000+(j+h)*M+i,
172
                                   Acomm );
173
                  HPL_dlaprnt( ib, 1, buf, i+1, j+h+1, NB, CMATNM );
174
               }
175
            }
176
            if( myrow == icurrow ) ii += ib;
177
            icurrow = MModAdd1( icurrow, nprow );
178
            (void) HPL_barrier( Acomm );
179
         }
180
         ii = 0; icurrow = IAROW;
181
      }
182
      if( mycol == icurcol ) jj += jb;
183
      icurcol = MModAdd1( icurcol, npcol );
184
      (void) HPL_barrier( Acomm );
185
   }
186
   if( ( myrow == 0 ) && ( mycol == 0 ) && ( buf ) ) free( buf );
187
/*
188
 * End of HPL_pdlaprnt
189
 */
190
}