Actual source code: ispai.c

  1: /*
  2:    3/99 Modified by Stephen Barnard to support SPAI version 3.0
  3: */

  5: /*
  6:       Provides an interface to the SPAI Sparse Approximate Inverse Preconditioner
  7:    Code written by Stephen Barnard.

  9:       Note: there is some BAD memory bleeding below!

 11:       This code needs work

 13:    1) get rid of all memory bleeding
 14:    2) fix PETSc/interface so that it gets if the matrix is symmetric from the matrix
 15:       rather than having the sp flag for PC_SPAI
 16:    3) fix to set the block size based on the matrix block size

 18: */
 19: #include <petsc/private/petscimpl.h>
 20: #include <petscpc.h>

 22: /*@
 23:   PCSPAISetEpsilon -- Set the tolerance for the `PCSPAI` preconditioner

 25:   Input Parameters:
 26: + pc       - the preconditioner
 27: - epsilon1 - the tolerance (default .4)

 29:   Level: intermediate

 31:   Note:
 32:   `espilon1` must be between 0 and 1. It controls the
 33:   quality of the approximation of M to the inverse of
 34:   A. Higher values of `epsilon1` lead to more work, more
 35:   fill, and usually better preconditioners. In many
 36:   cases the best choice of `epsilon1` is the one that
 37:   divides the total solution time equally between the
 38:   preconditioner and the solver.

 40: .seealso: [](ch_ksp), `PCSPAI`, `PCSetType()`
 41:   @*/
 42: PetscErrorCode PCSPAISetEpsilon(PC pc, PetscReal epsilon1)
 43: {
 44:   PetscFunctionBegin;
 45:   PetscTryMethod(pc, "PCSPAISetEpsilon_C", (PC, PetscReal), (pc, epsilon1));
 46:   PetscFunctionReturn(PETSC_SUCCESS);
 47: }

 49: /*@
 50:   PCSPAISetNBSteps - set maximum number of improvement steps per row in
 51:   the `PCSPAI` preconditioner

 53:   Input Parameters:
 54: + pc       - the preconditioner
 55: - nbsteps1 - number of steps (default 5)

 57:   Note:
 58:   `PCSPAI` constructs to approximation to every column of
 59:   the exact inverse of A in a series of improvement
 60:   steps. The quality of the approximation is determined
 61:   by epsilon. If an approximation achieving an accuracy
 62:   of epsilon is not obtained after `nbsteps1` steps, `PCSPAI` simply
 63:   uses the best approximation constructed so far.

 65:   Level: intermediate

 67: .seealso: [](ch_ksp), `PCSPAI`, `PCSetType()`, `PCSPAISetMaxNew()`
 68: @*/
 69: PetscErrorCode PCSPAISetNBSteps(PC pc, PetscInt nbsteps1)
 70: {
 71:   PetscFunctionBegin;
 72:   PetscTryMethod(pc, "PCSPAISetNBSteps_C", (PC, PetscInt), (pc, nbsteps1));
 73:   PetscFunctionReturn(PETSC_SUCCESS);
 74: }

 76: /* added 1/7/99 g.h. */
 77: /*@
 78:   PCSPAISetMax - set the size of various working buffers in the `PCSPAI` preconditioner

 80:   Input Parameters:
 81: + pc   - the preconditioner
 82: - max1 - size (default is 5000)

 84:   Level: intermediate

 86: .seealso: [](ch_ksp), `PCSPAI`, `PCSetType()`
 87: @*/
 88: PetscErrorCode PCSPAISetMax(PC pc, PetscInt max1)
 89: {
 90:   PetscFunctionBegin;
 91:   PetscTryMethod(pc, "PCSPAISetMax_C", (PC, PetscInt), (pc, max1));
 92:   PetscFunctionReturn(PETSC_SUCCESS);
 93: }

 95: /*@
 96:   PCSPAISetMaxNew - set maximum number of new nonzero candidates per step in the `PCSPAI` preconditioner

 98:   Input Parameters:
 99: + pc      - the preconditioner
100: - maxnew1 - maximum number (default 5)

102:   Level: intermediate

104: .seealso: [](ch_ksp), `PCSPAI`, `PCSetType()`, `PCSPAISetNBSteps()`
105: @*/
106: PetscErrorCode PCSPAISetMaxNew(PC pc, PetscInt maxnew1)
107: {
108:   PetscFunctionBegin;
109:   PetscTryMethod(pc, "PCSPAISetMaxNew_C", (PC, PetscInt), (pc, maxnew1));
110:   PetscFunctionReturn(PETSC_SUCCESS);
111: }

113: /*@
114:   PCSPAISetBlockSize - set the block size for the `PCSPAI` preconditioner

116:   Input Parameters:
117: + pc          - the preconditioner
118: - block_size1 - block size (default 1)

120:   Level: intermediate

122:   Notes:
123:   A block
124:   size of 1 treats A as a matrix of scalar elements. A
125:   block size of s > 1 treats A as a matrix of sxs
126:   blocks. A block size of 0 treats A as a matrix with
127:   variable sized blocks, which are determined by
128:   searching for dense square diagonal blocks in A.
129:   This can be very effective for finite-element
130:   matrices.

132:   SPAI will convert A to block form, use a block
133:   version of the preconditioner algorithm, and then
134:   convert the result back to scalar form.

136:   In many cases the a block-size parameter other than 1
137:   can lead to very significant improvement in
138:   performance.

140:   Developer Note:
141:   This preconditioner could use the matrix block size as the default block size to use

143: .seealso: [](ch_ksp), `PCSPAI`, `PCSetType()`
144: @*/
145: PetscErrorCode PCSPAISetBlockSize(PC pc, PetscInt block_size1)
146: {
147:   PetscFunctionBegin;
148:   PetscTryMethod(pc, "PCSPAISetBlockSize_C", (PC, PetscInt), (pc, block_size1));
149:   PetscFunctionReturn(PETSC_SUCCESS);
150: }

152: /*@
153:   PCSPAISetCacheSize - specify cache size in the `PCSPAI` preconditioner

155:   Input Parameters:
156: + pc         - the preconditioner
157: - cache_size - cache size {0,1,2,3,4,5} (default 5)

159:   Level: intermediate

161:   Note:
162:   `PCSPAI` uses a hash table to cache messages and avoid
163:   redundant communication. If suggest always using
164:   5. This parameter is irrelevant in the serial
165:   version.

167: .seealso: [](ch_ksp), `PCSPAI`, `PCSetType()`
168: @*/
169: PetscErrorCode PCSPAISetCacheSize(PC pc, PetscInt cache_size)
170: {
171:   PetscFunctionBegin;
172:   PetscTryMethod(pc, "PCSPAISetCacheSize_C", (PC, PetscInt), (pc, cache_size));
173:   PetscFunctionReturn(PETSC_SUCCESS);
174: }

176: /*@
177:   PCSPAISetVerbose - verbosity level for the `PCSPAI` preconditioner

179:   Input Parameters:
180: + pc      - the preconditioner
181: - verbose - level (default 1)

183:   Level: intermediate

185:   Note:
186:   Prints parameters, timings and matrix statistics

188: .seealso: [](ch_ksp), `PCSPAI`, `PCSetType()`
189: @*/
190: PetscErrorCode PCSPAISetVerbose(PC pc, PetscInt verbose)
191: {
192:   PetscFunctionBegin;
193:   PetscTryMethod(pc, "PCSPAISetVerbose_C", (PC, PetscInt), (pc, verbose));
194:   PetscFunctionReturn(PETSC_SUCCESS);
195: }

197: /*@
198:   PCSPAISetSp - specify a symmetric matrix sparsity pattern in the `PCSPAI` preconditioner

200:   Input Parameters:
201: + pc - the preconditioner
202: - sp - 0 or 1

204:   Level: intermediate

206:   Note:
207:   If A has a symmetric nonzero pattern use `sp` 1 to
208:   improve performance by eliminating some communication
209:   in the parallel version. Even if A does not have a
210:   symmetric nonzero pattern `sp` 1 may well lead to good
211:   results, but the code will not follow the published
212:   SPAI algorithm exactly.

214: .seealso: [](ch_ksp), `PCSPAI`, `PCSetType()`
215: @*/
216: PetscErrorCode PCSPAISetSp(PC pc, PetscInt sp)
217: {
218:   PetscFunctionBegin;
219:   PetscTryMethod(pc, "PCSPAISetSp_C", (PC, PetscInt), (pc, sp));
220:   PetscFunctionReturn(PETSC_SUCCESS);
221: }