PLASMA  2.4.5
PLASMA - Parallel Linear Algebra for Scalable Multi-core Architectures
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Groups
zgeqrs.c
Go to the documentation of this file.
1 
15 #include "common.h"
16 
17 /***************************************************************************/
67 int PLASMA_zgeqrs(int M, int N, int NRHS,
68  PLASMA_Complex64_t *A, int LDA,
70  PLASMA_Complex64_t *B, int LDB)
71 {
72  int NB, IB, IBNB, MT, NT;
73  int status;
75  PLASMA_sequence *sequence = NULL;
77  PLASMA_desc descA, descB, descT;
78 
79  plasma = plasma_context_self();
80  if (plasma == NULL) {
81  plasma_fatal_error("PLASMA_zgeqrs", "PLASMA not initialized");
83  }
84  /* Check input arguments */
85  if (M < 0) {
86  plasma_error("PLASMA_zgeqrs", "illegal value of M");
87  return -1;
88  }
89  if (N < 0 || N > M) {
90  plasma_error("PLASMA_zgeqrs", "illegal value of N");
91  return -2;
92  }
93  if (NRHS < 0) {
94  plasma_error("PLASMA_zgeqrs", "illegal value of N");
95  return -3;
96  }
97  if (LDA < max(1, M)) {
98  plasma_error("PLASMA_zgeqrs", "illegal value of LDA");
99  return -5;
100  }
101  if (LDB < max(1, max(1, M))) {
102  plasma_error("PLASMA_zgeqrs", "illegal value of LDB");
103  return -8;
104  }
105  /* Quick return */
106  if (min(M, min(N, NRHS)) == 0) {
107  return PLASMA_SUCCESS;
108  }
109 
110  /* Tune NB & IB depending on M, N & NRHS; Set NBNBSIZE */
111  status = plasma_tune(PLASMA_FUNC_ZGELS, M, N, NRHS);
112  if (status != PLASMA_SUCCESS) {
113  plasma_error("PLASMA_zgeqrs", "plasma_tune() failed");
114  return status;
115  }
116 
117  /* Set MT, NT & NTRHS */
118  NB = PLASMA_NB;
119  IB = PLASMA_IB;
120  IBNB = IB*NB;
121  MT = (M%NB==0) ? (M/NB) : (M/NB+1);
122  NT = (N%NB==0) ? (N/NB) : (N/NB+1);
123 
124  plasma_sequence_create(plasma, &sequence);
125 
126  if (plasma->householder == PLASMA_FLAT_HOUSEHOLDER) {
127  descT = plasma_desc_init(
129  IB, NB, IBNB,
130  MT*IB, NT*NB, 0, 0, MT*IB, NT*NB);
131  }
132  else {
133  /* Double the size of T to accomodate the tree reduction phase */
134  descT = plasma_desc_init(
136  IB, NB, IBNB,
137  MT*IB, 2*NT*NB, 0, 0, MT*IB, 2*NT*NB);
138  }
139  descT.mat = T;
140 
142  plasma_zooplap2tile( descA, A, NB, NB, LDA, N, 0, 0, M, N , plasma_desc_mat_free(&(descA)) );
143  plasma_zooplap2tile( descB, B, NB, NB, LDB, NRHS, 0, 0, M, NRHS, plasma_desc_mat_free(&(descA)); plasma_desc_mat_free(&(descB)));
144  } else {
145  plasma_ziplap2tile( descA, A, NB, NB, LDA, N, 0, 0, M, N );
146  plasma_ziplap2tile( descB, B, NB, NB, LDB, NRHS, 0, 0, M, NRHS);
147  }
148 
149  /* Call the tile interface */
150  PLASMA_zgeqrs_Tile_Async(&descA, &descT, &descB, sequence, &request);
151 
153  plasma_zooptile2lap( descA, A, NB, NB, LDA, N );
154  plasma_zooptile2lap( descB, B, NB, NB, LDB, NRHS );
156  plasma_desc_mat_free(&descA);
157  plasma_desc_mat_free(&descB);
158  } else {
159  plasma_ziptile2lap( descA, A, NB, NB, LDA, N );
160  plasma_ziptile2lap( descB, B, NB, NB, LDB, NRHS );
162  }
163 
164  status = sequence->status;
165  plasma_sequence_destroy(plasma, sequence);
166  return status;
167 }
168 
169 /***************************************************************************/
207 {
209  PLASMA_sequence *sequence = NULL;
211  int status;
212 
213  plasma = plasma_context_self();
214  if (plasma == NULL) {
215  plasma_fatal_error("PLASMA_zgeqrs_Tile", "PLASMA not initialized");
217  }
218  plasma_sequence_create(plasma, &sequence);
219  PLASMA_zgeqrs_Tile_Async(A, T, B, sequence, &request);
221  status = sequence->status;
222  plasma_sequence_destroy(plasma, sequence);
223  return status;
224 }
225 
226 /***************************************************************************/
256  PLASMA_sequence *sequence, PLASMA_request *request)
257 {
258  PLASMA_desc descA = *A;
259  PLASMA_desc descT = *T;
260  PLASMA_desc descB = *B;
262 
263  plasma = plasma_context_self();
264  if (plasma == NULL) {
265  plasma_fatal_error("PLASMA_zgeqrs_Tile", "PLASMA not initialized");
267  }
268  if (sequence == NULL) {
269  plasma_fatal_error("PLASMA_zgeqrs_Tile", "NULL sequence");
270  return PLASMA_ERR_UNALLOCATED;
271  }
272  if (request == NULL) {
273  plasma_fatal_error("PLASMA_zgeqrs_Tile", "NULL request");
274  return PLASMA_ERR_UNALLOCATED;
275  }
276  /* Check sequence status */
277  if (sequence->status == PLASMA_SUCCESS)
278  request->status = PLASMA_SUCCESS;
279  else
280  return plasma_request_fail(sequence, request, PLASMA_ERR_SEQUENCE_FLUSHED);
281 
282  /* Check descriptors for correctness */
283  if (plasma_desc_check(&descA) != PLASMA_SUCCESS) {
284  plasma_error("PLASMA_zgeqrs_Tile", "invalid first descriptor");
285  return plasma_request_fail(sequence, request, PLASMA_ERR_ILLEGAL_VALUE);
286  }
287  if (plasma_desc_check(&descT) != PLASMA_SUCCESS) {
288  plasma_error("PLASMA_zgeqrs_Tile", "invalid second descriptor");
289  return plasma_request_fail(sequence, request, PLASMA_ERR_ILLEGAL_VALUE);
290  }
291  if (plasma_desc_check(&descB) != PLASMA_SUCCESS) {
292  plasma_error("PLASMA_zgeqrs_Tile", "invalid third descriptor");
293  return plasma_request_fail(sequence, request, PLASMA_ERR_ILLEGAL_VALUE);
294  }
295  /* Check input arguments */
296  if (descA.nb != descA.mb || descB.nb != descB.mb) {
297  plasma_error("PLASMA_zgeqrs_Tile", "only square tiles supported");
298  return plasma_request_fail(sequence, request, PLASMA_ERR_ILLEGAL_VALUE);
299  }
300  /* Quick return */
301 /*
302  if (min(M, min(N, NRHS)) == 0) {
303  return PLASMA_SUCCESS;
304  }
305 */
306  if (plasma->householder == PLASMA_FLAT_HOUSEHOLDER) {
310  PLASMA_desc, descA,
311  PLASMA_desc, descB,
312  PLASMA_desc, descT,
313  PLASMA_sequence*, sequence,
314  PLASMA_request*, request);
315  }
316  else {
317  plasma_dynamic_call_8(plasma_pzunmqrrh,
320  PLASMA_desc, descA,
321  PLASMA_desc, descB,
322  PLASMA_desc, descT,
324  PLASMA_sequence*, sequence,
325  PLASMA_request*, request);
326  }
327 
333  PLASMA_Complex64_t, 1.0,
334  PLASMA_desc, plasma_desc_submatrix(descA, 0, 0, descA.n, descA.n),
335  PLASMA_desc, plasma_desc_submatrix(descB, 0, 0, descA.n, descB.n),
336  PLASMA_sequence*, sequence,
337  PLASMA_request*, request);
338 
339  return PLASMA_SUCCESS;
340 }