Skip to content
Snippets Groups Projects
lin_alg.h 3.92 KiB
Newer Older
Kurt A. O'Hearn's avatar
Kurt A. O'Hearn committed
/*----------------------------------------------------------------------
  SerialReax - Reax Force Field Simulator
Kurt A. O'Hearn's avatar
Kurt A. O'Hearn committed

Kurt A. O'Hearn's avatar
Kurt A. O'Hearn committed
  Copyright (2010) Purdue University
  Hasan Metin Aktulga, haktulga@cs.purdue.edu
  Joseph Fogarty, jcfogart@mail.usf.edu
  Sagar Pandit, pandit@usf.edu
  Ananth Y Grama, ayg@cs.purdue.edu
Kurt A. O'Hearn's avatar
Kurt A. O'Hearn committed

Kurt A. O'Hearn's avatar
Kurt A. O'Hearn committed
  This program is free software; you can redistribute it and/or
  modify it under the terms of the GNU General Public License as
Kurt A. O'Hearn's avatar
Kurt A. O'Hearn committed
  published by the Free Software Foundation; either version 2 of
Kurt A. O'Hearn's avatar
Kurt A. O'Hearn committed
  the License, or (at your option) any later version.
Kurt A. O'Hearn's avatar
Kurt A. O'Hearn committed

Kurt A. O'Hearn's avatar
Kurt A. O'Hearn committed
  This program is distributed in the hope that it will be useful,
  but WITHOUT ANY WARRANTY; without even the implied warranty of
Kurt A. O'Hearn's avatar
Kurt A. O'Hearn committed
  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
Kurt A. O'Hearn's avatar
Kurt A. O'Hearn committed
  See the GNU General Public License for more details:
  <http://www.gnu.org/licenses/>.
  ----------------------------------------------------------------------*/

#ifndef __LIN_ALG_H_
#define __LIN_ALG_H_
Kurt A. O'Hearn's avatar
Kurt A. O'Hearn committed

#include "mytypes.h"

typedef enum
{
    LOWER = 0,
    UPPER = 1,
} TRIANGULARITY;

void Sort_Matrix_Rows( sparse_matrix * const );

void setup_sparse_approx_inverse( const sparse_matrix * const, sparse_matrix **, 
        sparse_matrix **, sparse_matrix **, sparse_matrix **, const real );
int Estimate_LU_Fill( const sparse_matrix * const, const real * const );

void Calculate_Droptol( const sparse_matrix * const,
        real * const, const real );

#if defined(HAVE_SUPERLU_MT)
real SuperLU_Factorize( const sparse_matrix * const,
        sparse_matrix * const, sparse_matrix * const );
#endif

real diag_pre_comp( const sparse_matrix * const, real * const );

real ICHOLT( const sparse_matrix * const, const real * const,
        sparse_matrix * const, sparse_matrix * const );

#if defined(TESTING)
real ICHOL_PAR( const sparse_matrix * const, const unsigned int,
        sparse_matrix * const, sparse_matrix * const );
#endif

real ILU_PAR( const sparse_matrix * const, const unsigned int,
        sparse_matrix * const, sparse_matrix * const );

real ILUT_PAR( const sparse_matrix * const, const real *,
        const unsigned int, sparse_matrix * const, sparse_matrix * const );
#if defined(HAVE_LAPACKE) || defined(HAVE_LAPACKE_MKL)
real sparse_approx_inverse( const sparse_matrix * const, const sparse_matrix * const,
        sparse_matrix ** );
#endif

void Transpose( const sparse_matrix * const, sparse_matrix * const );
void Transpose_I( sparse_matrix * const );

void tri_solve( const sparse_matrix * const, const real * const,
        real * const, const int, const TRIANGULARITY );
void tri_solve_level_sched( static_storage *,
        const sparse_matrix * const,
        const real * const, real * const, const int,
        const TRIANGULARITY, int );
void jacobi_iter( const static_storage * const,
        const sparse_matrix * const, const real * const,
        const real * const, real * const, const TRIANGULARITY,
        const unsigned int );

void setup_graph_coloring( const control_params * const,
        const static_storage * const, const sparse_matrix * const,
int GMRES( const static_storage * const, const control_params * const,
        simulation_data * const, const sparse_matrix * const,
        const real * const, const real, real * const,
Kurt A. O'Hearn's avatar
Kurt A. O'Hearn committed

int GMRES_HouseHolder( const static_storage * const, const control_params * const,
        simulation_data * const, const sparse_matrix * const,
        const real * const, const real, real * const,
Kurt A. O'Hearn's avatar
Kurt A. O'Hearn committed

int CG( const static_storage * const, const control_params * const,
        simulation_data * const, const sparse_matrix * const, const real * const,
        const real, real * const, const int );
Kurt A. O'Hearn's avatar
Kurt A. O'Hearn committed

int SDM( const static_storage * const, const control_params * const,
        simulation_data * const, const sparse_matrix * const, const real * const, const real,
Kurt A. O'Hearn's avatar
Kurt A. O'Hearn committed

real condest( const sparse_matrix * const, const sparse_matrix * const );
Kurt A. O'Hearn's avatar
Kurt A. O'Hearn committed
#endif