SuperLU 6.0.0
Macros | Functions
fgmr.c File Reference

ARMS preconditioned FGMRES. More...

#include <stdio.h>
#include <stdlib.h>
#include <string.h>
#include <math.h>
#include "slu_ddefs.h"
Include dependency graph for fgmr.c:

Macros

#define epsmac   1.0e-16
 

Functions

double ddot_ (int *, double[], int *, double[], int *)
 
double dnrm2_ (int *, double[], int *)
 
void daxpy_ (int *, double *, double[], int *, double[], int *)
 
double dcopy_ (int *, double[], int *, double[], int *)
 
int fgmr (int n, void(*matvec)(double, double[], double, double[]), void(*psolve)(int, double[], double[]), double *rhs, double *sol, double tol, int im, int *itmax, FILE *fits)
 Simple version of the ARMS preconditioned FGMRES algorithm. More...
 

Macro Definition Documentation

◆ epsmac

#define epsmac   1.0e-16

Function Documentation

◆ daxpy_()

void daxpy_ ( int *  ,
double *  ,
double  [],
int *  ,
double  [],
int *   
)

◆ dcopy_()

double dcopy_ ( int *  ,
double  [],
int *  ,
double  [],
int *   
)

◆ ddot_()

double ddot_ ( int *  ,
double  [],
int *  ,
double  [],
int *   
)

◆ dnrm2_()

double dnrm2_ ( int *  ,
double  [],
int *   
)

◆ fgmr()

int fgmr ( int  n,
void(*)(double, double[], double, double[])  matvec,
void(*)(int, double[], double[])  psolve,
double *  rhs,
double *  sol,
double  tol,
int  im,
int *  itmax,
FILE *  fits 
)

Y. S. Dec. 2000. – Apr. 2008

internal work arrays: vv = work array of length [im+1][n] (used to store the Arnoldi basis) hh = work array of length [im][im+1] (Householder matrix) z = work array of length [im][n] to store preconditioned vectors

Parameters
[in]nDimension of vectors and matrices.
[in]matvecOperation for matrix-vector multiplication.
[in]psolve(right) preconditioning operation. Can be a NULL pointer (GMRES without preconditioner)
[in]rhsReal vector of length n containing the right hand side.
[in,out]solIn: Real vector of length n containing an initial guess to the solution on input. Out: Contains an approximate solution (upon successful return).
[in]tolTolerance for stopping iteration
[in]imKrylov subspace dimension
[in,out]itmaxIn: max number of iterations allowed. Out: number of steps required to converge.
[in]fitsIf NULL, no output. If not NULL, file handle to output "resid vs time and its".
Returns
Whether the algorithm finished successfully.
Here is the call graph for this function: