hypre/seq_mv/HYPRE_csr_matrix.c
falgout e3181f26b1 Added 64 bit feature using HYPRE_Int (see tracker [issue489] for details).
Changed MPI routines to hypre_MPI routines.
Added hypre_printf, etc. routines.
Added AUTOTEST tests to look for 'int' and 'MPI_' calls.
Added a new approach for the Fortran interface (not implemented everywhere yet).
2010-12-20 19:27:44 +00:00

107 lines
3.2 KiB
C

/*BHEADER**********************************************************************
* Copyright (c) 2008, Lawrence Livermore National Security, LLC.
* Produced at the Lawrence Livermore National Laboratory.
* This file is part of HYPRE. See file COPYRIGHT for details.
*
* HYPRE is free software; you can redistribute it and/or modify it under the
* terms of the GNU Lesser General Public License (as published by the Free
* Software Foundation) version 2.1 dated February 1999.
*
* $Revision$
***********************************************************************EHEADER*/
/******************************************************************************
*
* HYPRE_CSRMatrix interface
*
*****************************************************************************/
#include "headers.h"
/*--------------------------------------------------------------------------
* HYPRE_CSRMatrixCreate
*--------------------------------------------------------------------------*/
HYPRE_CSRMatrix
HYPRE_CSRMatrixCreate( HYPRE_Int num_rows,
HYPRE_Int num_cols,
HYPRE_Int *row_sizes )
{
hypre_CSRMatrix *matrix;
HYPRE_Int *matrix_i;
HYPRE_Int i;
matrix_i = hypre_CTAlloc(HYPRE_Int, num_rows + 1);
matrix_i[0] = 0;
for (i = 0; i < num_rows; i++)
{
matrix_i[i+1] = matrix_i[i] + row_sizes[i];
}
matrix = hypre_CSRMatrixCreate(num_rows, num_cols, matrix_i[num_rows]);
hypre_CSRMatrixI(matrix) = matrix_i;
return ( (HYPRE_CSRMatrix) matrix );
}
/*--------------------------------------------------------------------------
* HYPRE_CSRMatrixDestroy
*--------------------------------------------------------------------------*/
HYPRE_Int
HYPRE_CSRMatrixDestroy( HYPRE_CSRMatrix matrix )
{
return( hypre_CSRMatrixDestroy( (hypre_CSRMatrix *) matrix ) );
}
/*--------------------------------------------------------------------------
* HYPRE_CSRMatrixInitialize
*--------------------------------------------------------------------------*/
HYPRE_Int
HYPRE_CSRMatrixInitialize( HYPRE_CSRMatrix matrix )
{
return ( hypre_CSRMatrixInitialize( (hypre_CSRMatrix *) matrix ) );
}
/*--------------------------------------------------------------------------
* HYPRE_CSRMatrixRead
*--------------------------------------------------------------------------*/
HYPRE_CSRMatrix
HYPRE_CSRMatrixRead( char *file_name )
{
return ( (HYPRE_CSRMatrix) hypre_CSRMatrixRead( file_name ) );
}
/*--------------------------------------------------------------------------
* HYPRE_CSRMatrixPrint
*--------------------------------------------------------------------------*/
void
HYPRE_CSRMatrixPrint( HYPRE_CSRMatrix matrix,
char *file_name )
{
hypre_CSRMatrixPrint( (hypre_CSRMatrix *) matrix,
file_name );
}
/*--------------------------------------------------------------------------
* HYPRE_CSRMatrixGetNumRows
*--------------------------------------------------------------------------*/
HYPRE_Int
HYPRE_CSRMatrixGetNumRows( HYPRE_CSRMatrix matrix, HYPRE_Int *num_rows )
{
hypre_CSRMatrix *csr_matrix = (hypre_CSRMatrix *) matrix;
*num_rows = hypre_CSRMatrixNumRows( csr_matrix );
return 0;
}