StarNEig User's Guide  v0.1.7
A task-based library for solving dense nonsymmetric eigenvalue problems
sep_dm_full_chain.c
#include "validate.h"
#include <stdlib.h>
#include <stdio.h>
#include <time.h>
#include <mpi.h>
// a predicate function that selects all eigenvalues that have positive a real
// part
static int predicate(double real, double imag, void *arg)
{
if (0.0 < real)
return 1;
return 0;
}
int main(int argc, char **argv)
{
const int n = 3000; // matrix dimension
const int root = 0; // root rank
// initialize MPI
int thread_support;
MPI_Init_thread(
&argc, (char ***)&argv, MPI_THREAD_MULTIPLE, &thread_support);
int world_rank;
MPI_Comm_rank(MPI_COMM_WORLD, &world_rank);
// the root node initializes the matrices locally
int ldA = 0, ldQ = 0, ldC = 0;
double *A = NULL, *Q = NULL, *C = NULL;
if (world_rank == root) {
srand((unsigned) time(NULL));
// generate a full random matrix A and a copy C
ldA = ((n/8)+1)*8; ldC = ((n/8)+1)*8;
A = malloc(n*ldA*sizeof(double));
C = malloc(n*ldC*sizeof(double));
for (int j = 0; j < n; j++)
for (int i = 0; i < n; i++)
A[j*ldA+i] = C[j*ldC+i] = 2.0*rand()/RAND_MAX - 1.0;
// generate an identity matrix Q
ldQ = ((n/8)+1)*8;
Q = malloc(n*ldA*sizeof(double));
for (int j = 0; j < n; j++)
for (int i = 0; i < n; i++)
Q[j*ldQ+i] = i == j ? 1.0 : 0.0;
}
// allocate space for the eigenvalues and the eigenvalue selection vector
double *real = malloc(n*sizeof(double));
double *imag = malloc(n*sizeof(double));
int *select = malloc(n*sizeof(int));
// Initialize the StarNEig library using all available CPU cores and
// GPUs. The STARNEIG_HINT_DM flag indicates that the library should
// initialize itself for distributed memory computations.
// create a two-dimensional block cyclic distribution with row-major
// ordering
// Convert the local matrix A to a distributed matrix lA that is owned by
// the root node. This is done in-place, i.e., the matrices A and lA point
// to the same data.
n, n, STARNEIG_REAL_DOUBLE, root, A, ldA);
// create a distributed matrix dA using default data distribution and
// distributed block size
// copy the local matrix lA to the distributed matrix dA (scatter)
// scatter the matrix Q
n, n, STARNEIG_REAL_DOUBLE, root, Q, ldQ);
// reduce the full matrix dA to upper Hessenberg form
printf("Hessenberg reduction...\n");
// reduce the upper Hessenberg matrix dA to Schur form
printf("Schur reduction...\n");
starneig_SEP_DM_Schur(dA, dQ, real, imag);
// select eigenvalues that have positive a real part
int num_selected;
starneig_SEP_DM_Select(dA, &predicate, NULL, select, &num_selected);
printf("Selected %d eigenvalues out of %d.\n", num_selected, n);
// reorder the selected eigenvalues to the upper left corner of the matrix
// dA
printf("Reordering...\n");
starneig_SEP_DM_ReorderSchur(select, dA, dQ, real, imag);
// copy the distributed matrix dA back to the local matrix lA (gather)
// free the distributed matrix lA (matrix A is not freed)
// free the distributed matrix dA (all local resources are freed)
// gather the matrix Q
// free the data distribution
// de-initialize the StarNEig library
// de-initialize MPI
MPI_Finalize();
if (world_rank == root) {
// check residual || Q A Q^T - C ||_F / || C ||_F
check_residual(n, ldQ, ldA, ldQ, ldC, Q, A, Q, C);
// check residual || Q Q^T - I ||_F / || I ||_F
check_orthogonality(n, ldQ, Q);
}
// cleanup
free(A);
free(C);
free(Q);
free(real);
free(imag);
free(select);
return 0;
}
STARNEIG_ORDER_ROW_MAJOR
@ STARNEIG_ORDER_ROW_MAJOR
Row-major natural ordering.
Definition: distr_matrix.h:80
starneig_SEP_DM_ReorderSchur
starneig_error_t starneig_SEP_DM_ReorderSchur(int selected[], starneig_distr_matrix_t S, starneig_distr_matrix_t Q, double real[], double imag[])
Reorders selected eigenvalues to the top left corner of a Schur decomposition.
starneig.h
This file includes most StarNEig header files.
starneig_distr_matrix_create_local
starneig_distr_matrix_t starneig_distr_matrix_create_local(int rows, int cols, starneig_datatype_t type, int owner, double *A, int ldA)
Creates a single-owner distributed matrix from a local matrix.
starneig_distr_init_mesh
starneig_distr_t starneig_distr_init_mesh(int rows, int cols, starneig_distr_order_t order)
Creates a two-dimensional block cyclic data distribution.
starneig_distr_matrix_destroy
void starneig_distr_matrix_destroy(starneig_distr_matrix_t matrix)
Destroys a distributed matrix.
starneig_distr_destroy
void starneig_distr_destroy(starneig_distr_t distr)
Destroys a data distribution.
STARNEIG_USE_ALL
#define STARNEIG_USE_ALL
Use all resources.
Definition: node.h:72
starneig_SEP_DM_Hessenberg
starneig_error_t starneig_SEP_DM_Hessenberg(starneig_distr_matrix_t A, starneig_distr_matrix_t Q)
Computes a Hessenberg decomposition of a general matrix.
starneig_distr_matrix_copy
void starneig_distr_matrix_copy(starneig_distr_matrix_t source, starneig_distr_matrix_t dest)
Copies the contents of a distributed matrix to a second distributed matrix.
STARNEIG_HINT_DM
#define STARNEIG_HINT_DM
Distributed memory mode.
Definition: node.h:102
starneig_node_init
void starneig_node_init(int cores, int gpus, starneig_flag_t flags)
Initializes the intra-node execution environment.
starneig_SEP_DM_Schur
starneig_error_t starneig_SEP_DM_Schur(starneig_distr_matrix_t H, starneig_distr_matrix_t Q, double real[], double imag[])
Computes a Schur decomposition given a Hessenberg decomposition.
starneig_distr_matrix_t
struct starneig_distr_matrix * starneig_distr_matrix_t
Distributed matrix.
Definition: distr_matrix.h:189
STARNEIG_REAL_DOUBLE
@ STARNEIG_REAL_DOUBLE
Double precision real numbers.
Definition: distr_matrix.h:195
starneig_distr_matrix_create
starneig_distr_matrix_t starneig_distr_matrix_create(int rows, int cols, int row_blksz, int col_blksz, starneig_datatype_t type, starneig_distr_t distr)
Creates a distributed matrix with uninitialized matrix elements.
starneig_SEP_DM_Select
starneig_error_t starneig_SEP_DM_Select(starneig_distr_matrix_t S, int(*predicate)(double real, double imag, void *arg), void *arg, int selected[], int *num_selected)
Generates a selection array for a Schur matrix using a user-supplied predicate function.
starneig_node_finalize
void starneig_node_finalize()
Deallocates resources associated with the intra-node configuration.
starneig_distr_t
struct starneig_distr * starneig_distr_t
Data distribution.
Definition: distr_matrix.h:73