How to use LAPACK on maya

Introduction

The Linear Algebra PACKage (LAPACK) library is useful for matrix algebra operations like singular value decomposition, solving systems of linear equations, and computing eigenvalues. You may want to consider using it instead of writing your own routines. For lower level operations like matrix multiplication, see BLAS. In this tutorial, we will compile and run a program that uses LAPACK to compute SVDs. Before you begin, make sure to read the tutorial for compiling C programs.

Where to find documentation

On this page, we’ll show how to use LAPACK functionality within AMD Core Math Library (ACML). We’ll assume that you’re using the default module settings but make sure that you have the “acml/intel/64” module loaded.

For more information about the LAPACK routines in ACML, see

  • AMD’s ACML website
  • LAPACK at Netlib has a lot of information about the library. One useful section contains FORTRAN code for the procedures
  • The manual pages on maya (e.g. “man dgesvd”) also show the FORTRAN function interfaces
  • LAPACK User’s Guide
  • The ACML header file on maya shows the C interface. This file is located in: /cm/shared/apps/acml/5.3.1/ifort64/include/acml.h. You can browse through this file and see all the available functions and their arguments & return values. Note that there are two versions of most functions listed – one with an underscore (“_”) at the end, and one without an underscore. Both are versions of the same functional, but the underscore version is compatible with FORTRAN conventions, and the no-underscore version is a bit easier to use.

Example

In this example we will create a matrix, then compute its singular value decomposition (SVD) using the LAPACK dgesvd function. The SVD of A (m x n) will consist of U (m x m), S (p x 1, where p = min(m,n)), and VT (n x n) Our matrices will contain doubles, and be stored in column-major order. Here is the code

#include <stdio.h>
#include <acml.h>

#define MATRIX_IDX(n, i, j) j*n + i
#define MATRIX_ELEMENT(A, m, n, i, j) A[ MATRIX_IDX(m, i, j) ]

void init_matrix(double* A, int m, int n)
{
   for (int j = 0; j < n; j++)
   {
      for (int i = 0; i < m; i++)
      {
         MATRIX_ELEMENT(A, m, n, i, j) = 1.0 / (i + j + 1);
      }
   }
}

void print_matrix(const double* A, int m, int n)
{
   for (int i = 0; i < m; i++)
   {
      for (int j = 0; j < n; j++)
      {
          printf("%8.4f", MATRIX_ELEMENT(A, m, n, i, j));
      }
      printf("\n");
   }
}

int main(int argc, char** argv)
{
   int m = 3;
   int n = 4;

   int p = (m < n ? m : n);

   double A[m * n];
   double U[m * m];
   double VT[n * n];
   double S[p * 1];

   int info;

   init_matrix(A, m, n);

   printf("Matrix A (%d x %d) is:\n", m, n);
   print_matrix(A, m, n);

   dgesvd('A', 'A', m, n, A, m, S, U, m, VT, n, &info);
   if (info != 0)
   {
      fprintf(stderr, "Warning: dgesvd returned with a non-zero status (info = %d)\n", info);
   }

   printf("\nMatrix U (%d x %d) is:\n", m, m);
   print_matrix(U, m, m);

   printf("\nVector S (%d x %d) is:\n", p, 1);
   print_matrix(S, p, 1);

   printf("\nMatrix VT (%d x %d) is:\n", n, n);
   print_matrix(VT, n, n);

   return 0;
}


Download: ../code/lapack-svd/acml/main.c

The important line to notice, where the actual SVD computation is happening, is this

dgesvd('A', 'A', m, n, A, m, S, U, m, VT, n, &info);

where the ‘A’ arguments are flags saying that we want the entire U and VT matrices filled in by the procedure. For more information, see the documentation at Netlib, or try “man dgesvd” at the command line. Here is the Makefile

EXECUTABLE := matrix_svd
OBJS := main.o

CC := mpiicc
CFLAGS := -O3 -std=c99

INCLUDES :=
LIBLOCS :=
LDFLAGS := -lm -lacml

%.o: %.c %.h
    $(CC) $(CFLAGS) $(DEFS) $(INCLUDES) -c $< -o $@

$(EXECUTABLE): $(OBJS)
    $(CC) $(CFLAGS) $(DEFS) $(INCLUDES) $(OBJS) -o $@ $(LIBLOCS) $(LDFLAGS)

clean:
    -rm -f *.o $(EXECUTABLE)


Download: ../code/lapack-svd/acml/Makefile

The important part is that we need to link to two more libraries: libacml and libpgftnrtl. This is accomplished by adding them to the LDFLAGS variable. Compiling the code should look something like this

[araim1@maya-usr1 acml]$ make
mpiicc -O3 -std=c99   -c -o main.o main.c
mpiicc -O3 -std=c99   main.o -o matrix_svd  -lm -lacml
[araim1@maya-usr1 acml]$ ls
main.c  main.o  Makefile  matrix_svd
[araim1@maya-usr1 acml]$ 

Then running the code should look this

[araim1@maya-usr1 acml]$ ./matrix_svd
Matrix A (3 x 4) is:
  1.0000  0.5000  0.3333  0.2500
  0.5000  0.3333  0.2500  0.2000
  0.3333  0.2500  0.2000  0.1667

Matrix U (3 x 3) is:
 -0.8199  0.5563  0.1349
 -0.4662 -0.5123 -0.7213
 -0.3322 -0.6543  0.6794

Vector S (3 x 1) is:
  1.4519
  0.1433
  0.0042

Matrix VT (4 x 4) is:
 -0.8015 -0.4466 -0.3143 -0.2435
  0.5729 -0.3919 -0.5127 -0.5053
  0.1692 -0.7398  0.1245  0.6392
 -0.0263  0.3157 -0.7892  0.5261
[araim1@maya-usr1 acml]$ 

Running the code through the batch system should require no special steps. A standard SLURM script like the one below should be sufficient.

#!/bin/bash
#SBATCH --job-name=matrix_svd
#SBATCH --output=slurm.out
#SBATCH --error=slurm.err
#SBATCH --partition=develop
#SBATCH --nodes=1
#SBATCH --ntasks-per-node=1

./matrix_svd

Download: ../code/lapack-svd/acml/run.slurm