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
example_sormqr.c File Reference

Example for solving overdetermined linear systems with PLASMA_sormqr. More...

#include <stdlib.h>
#include <stdio.h>
#include <string.h>
#include <math.h>
#include <plasma.h>
#include <cblas.h>
#include <lapacke.h>
#include <core_blas.h>
Include dependency graph for example_sormqr.c:

Go to the source code of this file.

Macros

#define max(a, b)   ((a) > (b) ? (a) : (b))
#define min(a, b)   ((a) < (b) ? (a) : (b))

Functions

int check_solution (int, int, int, float *, int, float *, float *, int)
int main ()

Variables

int IONE = 1
int ISEED [4] = {0,0,0,1}

Detailed Description

Example for solving overdetermined linear systems with PLASMA_sormqr.

PLASMA testing routines PLASMA is a software package provided by Univ. of Tennessee, Univ. of California Berkeley and Univ. of Colorado Denver

Version:
2.4.5
Author:
Bilel Hadri
Date:
2010-11-15 s Tue Nov 22 14:35:54 2011

Definition in file example_sormqr.c.


Macro Definition Documentation

#define max (   a,
 
)    ((a) > (b) ? (a) : (b))

Definition at line 28 of file example_sormqr.c.

#define min (   a,
 
)    ((a) < (b) ? (a) : (b))

Definition at line 31 of file example_sormqr.c.


Function Documentation

int check_solution ( int  M,
int  N,
int  NRHS,
float *  A1,
int  LDA,
float *  B1,
float *  B2,
int  LDB 
)
int main ( )

Definition at line 39 of file example_sormqr.c.

References check_solution(), IONE, ISEED, PLASMA_Alloc_Workspace_sgeqrf(), PLASMA_Finalize(), PLASMA_Init(), PLASMA_sgeqrf(), PLASMA_sormqr(), PLASMA_strsm(), PlasmaLeft, PlasmaNonUnit, PlasmaNoTrans, PlasmaTrans, PlasmaUpper, and T.

{
int cores = 2;
int M = 15;
int N = 10;
int LDA = 15;
int NRHS = 5;
int LDB = 15;
int info;
int info_solution;
int i,j;
int LDAxN = LDA*N;
int LDBxNRHS = LDB*NRHS;
float *A1 = (float *)malloc(LDA*N*sizeof(float));
float *A2 = (float *)malloc(LDA*N*sizeof(float));
float *B1 = (float *)malloc(LDB*NRHS*sizeof(float));
float *B2 = (float *)malloc(LDB*NRHS*sizeof(float));
float *T;
/* Check if unable to allocate memory */
if ((!A1)||(!A2)||(!B1)||(!B2)){
printf("Out of Memory \n ");
return EXIT_SUCCESS;
}
/* Plasma Initialization */
PLASMA_Init(cores);
printf("-- PLASMA is initialized to run on %d cores. \n",cores);
/* Allocate T */
/* Initialize A1 and A2 */
LAPACKE_slarnv_work(IONE, ISEED, LDAxN, A1);
for (i = 0; i < M; i++)
for (j = 0; j < N; j++)
A2[LDA*j+i] = A1[LDA*j+i] ;
/* Initialize B1 and B2 */
LAPACKE_slarnv_work(IONE, ISEED, LDBxNRHS, B1);
for (i = 0; i < M; i++)
for (j = 0; j < NRHS; j++)
B2[LDB*j+i] = B1[LDB*j+i] ;
/* Factorization QR of the matrix A2 */
info = PLASMA_sgeqrf(M, N, A2, LDA, T);
/* Solve the problem */
info = PLASMA_sormqr(PlasmaLeft, PlasmaTrans, M, NRHS, N, A2, LDA, T, B2, LDB);
N, NRHS, (float)1.0, A2, LDA, B2, LDB);
/* Check the solution */
info_solution = check_solution(M, N, NRHS, A1, LDA, B1, B2, LDB);
if ((info_solution != 0)|(info != 0))
printf("-- Error in SORMQR example ! \n");
else
printf("-- Run of SORMQR example successful ! \n");
free(A1); free(A2); free(B1); free(B2); free(T);
return EXIT_SUCCESS;
}

Here is the call graph for this function:


Variable Documentation

int IONE = 1

Definition at line 36 of file example_sormqr.c.

int ISEED[4] = {0,0,0,1}

Definition at line 37 of file example_sormqr.c.