Developer Reference

Intel® oneAPI Math Kernel Library LAPACK Examples

ID 766877
Date 10/31/2024
Public
Document Table of Contents

LAPACKE_cheevr Example Program in C for Column Major Data Layout

/*******************************************************************************
*  Copyright (C) 2009-2015 Intel Corporation. All Rights Reserved.
*  The information and material ("Material") provided below is owned by Intel
*  Corporation or its suppliers or licensors, and title to such Material remains
*  with Intel Corporation or its suppliers or licensors. The Material contains
*  proprietary information of Intel or its suppliers and licensors. The Material
*  is protected by worldwide copyright laws and treaty provisions. No part of
*  the Material may be copied, reproduced, published, uploaded, posted,
*  transmitted, or distributed in any way without Intel's prior express written
*  permission. No license under any patent, copyright or other intellectual
*  property rights in the Material is granted to or conferred upon you, either
*  expressly, by implication, inducement, estoppel or otherwise. Any license
*  under such intellectual property rights must be express and approved by Intel
*  in writing.
*
********************************************************************************
*/
/*
   LAPACKE_cheevr Example.
   =======================

   Program computes eigenvalues specified by a selected range of values
   and corresponding eigenvectors of a complex Hermitian matrix A using the
   Relatively Robust Representations, where A is:

   ( -2.16,  0.00) ( -0.16, -4.86) ( -7.23, -9.38) ( -0.04,  6.86)
   ( -0.16,  4.86) (  7.45,  0.00) (  4.39,  6.29) ( -8.11, -4.41)
   ( -7.23,  9.38) (  4.39, -6.29) ( -9.03,  0.00) ( -6.89, -7.66)
   ( -0.04, -6.86) ( -8.11,  4.41) ( -6.89,  7.66) (  7.76,  0.00)

   Description.
   ============

   The routine computes selected eigenvalues and, optionally, eigenvectors of
   an n-by-n complex Hermitian matrix A. The eigenvector v(j) of A satisfies

   A*v(j) = lambda(j)*v(j)

   where lambda(j) is its eigenvalue. The computed eigenvectors are
   orthonormal.
   Eigenvalues and eigenvectors can be selected by specifying either a range
   of values or a range of indices for the desired eigenvalues.

   Example Program Results.
   ========================

 LAPACKE_cheevr (column-major, high-level) Example Program Results

 The total number of eigenvalues found: 2

 Selected eigenvalues
  -4.18   3.57

 Selected eigenvectors (stored columnwise)
 (  0.68,  0.00) (  0.38,  0.00)
 (  0.03,  0.18) (  0.54, -0.57)
 ( -0.03,  0.21) ( -0.40,  0.04)
 (  0.20,  0.64) ( -0.14, -0.26)
*/
#include <stdlib.h>
#include <stdio.h>
#include "mkl_lapacke.h"

/* Auxiliary routines prototypes */
extern void print_matrix( char* desc, MKL_INT m, MKL_INT n, MKL_Complex8* a, MKL_INT lda );
extern void print_rmatrix( char* desc, MKL_INT m, MKL_INT n, float* a, MKL_INT lda );

/* Parameters */
#define N 4
#define LDA N
#define LDZ N

/* Main program */
int main() {
        /* Locals */
        MKL_INT n = N, lda = LDA, ldz = LDZ, il, iu, m, info;
        float abstol, vl, vu;
        /* Local arrays */
        MKL_INT isuppz[N];
        float w[N];
        MKL_Complex8 z[LDZ*N];
        MKL_Complex8 a[LDA*N] = {
           {-2.16f,  0.00f}, {-0.16f,  4.86f}, {-7.23f,  9.38f}, {-0.04f, -6.86f},
           { 0.00f,  0.00f}, { 7.45f,  0.00f}, { 4.39f, -6.29f}, {-8.11f,  4.41f},
           { 0.00f,  0.00f}, { 0.00f,  0.00f}, {-9.03f,  0.00f}, {-6.89f,  7.66f},
           { 0.00f,  0.00f}, { 0.00f,  0.00f}, { 0.00f,  0.00f}, { 7.76f,  0.00f}
        };
        /* Executable statements */
        printf( "LAPACKE_cheevr (column-major, high-level) Example Program Results\n" );
        /* Negative abstol means using the default value */
        abstol = -1.0;
        /* Set VL, VU to compute eigenvalues in half-open (VL,VU] interval */
        vl = -5.0;
        vu = 5.0;
        /* Solve eigenproblem */
        info = LAPACKE_cheevr( LAPACK_COL_MAJOR, 'V', 'V', 'L', n, a, lda,
                        vl, vu, il, iu, abstol, &m, w, z, ldz, isuppz );
        /* Check for convergence */
        if( info > 0 ) {
                printf( "The algorithm failed to compute eigenvalues.\n" );
                exit( 1 );
        }
        /* Print the number of eigenvalues found */
        printf( "\n The total number of eigenvalues found:%2i\n", m );
        /* Print eigenvalues */
        print_rmatrix( "Selected eigenvalues", 1, m, w, 1 );
        /* Print eigenvectors */
        print_matrix( "Selected eigenvectors (stored columnwise)", n, m, z, ldz );
        exit( 0 );
} /* End of LAPACKE_cheevr Example */

/* Auxiliary routine: printing a matrix */
void print_matrix( char* desc, MKL_INT m, MKL_INT n, MKL_Complex8* a, MKL_INT lda ) {
        MKL_INT i, j;
        printf( "\n %s\n", desc );
        for( i = 0; i < m; i++ ) {
                for( j = 0; j < n; j++ )
                        printf( " (%6.2f,%6.2f)", a[i+j*lda].real, a[i+j*lda].imag );
                printf( "\n" );
        }
}

/* Auxiliary routine: printing a real matrix */
void print_rmatrix( char* desc, MKL_INT m, MKL_INT n, float* a, MKL_INT lda ) {
        MKL_INT i, j;
        printf( "\n %s\n", desc );
        for( i = 0; i < m; i++ ) {
                for( j = 0; j < n; j++ ) printf( " %6.2f", a[i+j*lda] );
                printf( "\n" );
        }
}