1 /****************************************************************************/ 2 /* mlsolv.h */ 3 /****************************************************************************/ 4 /* */ 5 /* Multi-Level SOLVers */ 6 /* */ 7 /* Copyright (C) 1992-1995 Tomas Skalicky. All rights reserved. */ 8 /* */ 9 /****************************************************************************/ 10 /* */ 11 /* ANY USE OF THIS CODE CONSTITUTES ACCEPTANCE OF THE TERMS */ 12 /* OF THE COPYRIGHT NOTICE (SEE FILE COPYRGHT.H) */ 13 /* */ 14 /****************************************************************************/ 15 16 #ifndef MLSOLV_H 17 #define MLSOLV_H 18 19 #include "laspack/vector.h" 20 #include "laspack/matrix.h" 21 #include "laspack/qmatrix.h" 22 #include "laspack/itersolv.h" 23 #include "laspack/copyrght.h" 24 25 Vector *MGStep(int NoLevels, QMatrix *A, Vector *x, Vector *b, 26 Matrix *R, Matrix *P, int Level, int Gamma, 27 IterProcType SmoothProc, int Nu1, int Nu2, 28 PrecondProcType PrecondProc, double Omega, 29 IterProcType SolvProc, int NuC, 30 PrecondProcType PrecondProcC, double OmegaC); 31 Vector *MGIter(int NoLevels, QMatrix *A, Vector *x, Vector *b, 32 Matrix *R, Matrix *P, int MaxIter, int Gamma, 33 IterProcType SmoothProc, int Nu1, int Nu2, 34 PrecondProcType PrecondProc, double Omega, 35 IterProcType SolvProc, int NuC, 36 PrecondProcType PrecondProcC, double OmegaC); 37 Vector *NestedMGIter(int NoLevels, QMatrix *A, Vector *x, Vector *b, 38 Matrix *R, Matrix *P, int Gamma, 39 IterProcType SmoothProc, int Nu1, int Nu2, 40 PrecondProcType PrecondProc, double Omega, 41 IterProcType SolvProc, int NuC, 42 PrecondProcType PrecondProcC, double OmegaC); 43 Vector *MGPCGIter(int NoLevels, QMatrix *A, Vector *x, Vector *b, 44 Matrix *R, Matrix *P, int MaxIter, int NoMGIter, int Gamma, 45 IterProcType SmoothProc, int Nu1, int Nu2, 46 PrecondProcType PrecondProc, double Omega, 47 IterProcType SolvProc, int NuC, 48 PrecondProcType PrecondProcC, double OmegaC); 49 Vector *BPXPrecond(int NoLevels, QMatrix *A, Vector *y, Vector *c, 50 Matrix *R, Matrix *P, int Level, 51 IterProcType SmoothProc, int Nu, 52 PrecondProcType PrecondProc, double Omega, 53 IterProcType SmoothProcC, int NuC, 54 PrecondProcType PrecondProcC, double OmegaC); 55 Vector *BPXPCGIter(int NoLevels, QMatrix *A, Vector *x, Vector *b, 56 Matrix *R, Matrix *P, int MaxIter, 57 IterProcType SmoothProc, int Nu, 58 PrecondProcType PrecondProc, double Omega, 59 IterProcType SmoothProcC, int NuC, 60 PrecondProcType PrecondProcC, double OmegaC); 61 62 #endif /* MLSOLV_H */ 63