Actual source code: mfnsolve.c

  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.
  7:    SLEPc is distributed under a 2-clause BSD license (see LICENSE).
  8:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  9: */
 10: /*
 11:    MFN routines related to the solution process
 12: */

 14: #include <slepc/private/mfnimpl.h>

 16: static PetscErrorCode MFNSolve_Private(MFN mfn,Vec b,Vec x)
 17: {
 18:   PetscFunctionBegin;
 19:   PetscCall(VecSetErrorIfLocked(x,3));

 21:   /* call setup */
 22:   PetscCall(MFNSetUp(mfn));
 23:   mfn->its = 0;

 25:   PetscCall(MFNViewFromOptions(mfn,NULL,"-mfn_view_pre"));

 27:   /* check nonzero right-hand side */
 28:   PetscCall(VecNorm(b,NORM_2,&mfn->bnorm));
 29:   PetscCheck(mfn->bnorm,PetscObjectComm((PetscObject)mfn),PETSC_ERR_ARG_WRONG,"Cannot pass a zero b vector to MFNSolve()");

 31:   /* call solver */
 32:   PetscCall(PetscLogEventBegin(MFN_Solve,mfn,b,x,0));
 33:   if (b!=x) PetscCall(VecLockReadPush(b));
 34:   PetscUseTypeMethod(mfn,solve,b,x);
 35:   if (b!=x) PetscCall(VecLockReadPop(b));
 36:   PetscCall(PetscLogEventEnd(MFN_Solve,mfn,b,x,0));

 38:   PetscCheck(mfn->reason,PetscObjectComm((PetscObject)mfn),PETSC_ERR_PLIB,"Internal error, solver returned without setting converged reason");

 40:   PetscCheck(!mfn->errorifnotconverged || mfn->reason>=0,PetscObjectComm((PetscObject)mfn),PETSC_ERR_NOT_CONVERGED,"MFNSolve has not converged");

 42:   /* various viewers */
 43:   PetscCall(MFNViewFromOptions(mfn,NULL,"-mfn_view"));
 44:   PetscCall(MFNConvergedReasonViewFromOptions(mfn));
 45:   PetscCall(MatViewFromOptions(mfn->A,(PetscObject)mfn,"-mfn_view_mat"));
 46:   PetscCall(VecViewFromOptions(b,(PetscObject)mfn,"-mfn_view_rhs"));
 47:   PetscCall(VecViewFromOptions(x,(PetscObject)mfn,"-mfn_view_solution"));
 48:   PetscFunctionReturn(PETSC_SUCCESS);
 49: }

 51: /*@
 52:    MFNSolve - Solves the matrix function problem. Given a vector b, the
 53:    vector x = f(A)*b is returned.

 55:    Collective

 57:    Input Parameters:
 58: +  mfn - the matrix function solver context
 59: -  b   - the right hand side vector

 61:    Output Parameter:
 62: .  x   - the solution (this may be the same vector as b, then b will be
 63:          overwritten with the answer)

 65:    Options Database Keys:
 66: +  -mfn_view - print information about the solver used
 67: .  -mfn_view_mat binary - save the matrix to the default binary viewer
 68: .  -mfn_view_rhs binary - save right hand side vector to the default binary viewer
 69: .  -mfn_view_solution binary - save computed solution vector to the default binary viewer
 70: -  -mfn_converged_reason - print reason for convergence, and number of iterations

 72:    Notes:
 73:    The matrix A is specified with MFNSetOperator().
 74:    The function f is specified with MFNSetFN().

 76:    Level: beginner

 78: .seealso: [](ch:mfn), `MFNCreate()`, `MFNSetUp()`, `MFNDestroy()`, `MFNSetTolerances()`, `MFNSetOperator()`, `MFNSetFN()`
 79: @*/
 80: PetscErrorCode MFNSolve(MFN mfn,Vec b,Vec x)
 81: {
 82:   PetscFunctionBegin;
 85:   PetscCheckSameComm(mfn,1,b,2);
 87:   if (b!=x) PetscCheckSameComm(mfn,1,x,3);
 88:   mfn->transpose_solve = PETSC_FALSE;
 89:   PetscCall(MFNSolve_Private(mfn,b,x));
 90:   PetscFunctionReturn(PETSC_SUCCESS);
 91: }

 93: /*@
 94:    MFNSolveTranspose - Solves the transpose matrix function problem. Given a vector b,
 95:    the vector x = f(A^T)*b is returned.

 97:    Collective

 99:    Input Parameters:
100: +  mfn - the matrix function solver context
101: -  b   - the right hand side vector

103:    Output Parameter:
104: .  x   - the solution (this may be the same vector as b, then b will be
105:          overwritten with the answer)

107:    Note:
108:    See available options at MFNSolve().

110:    Level: beginner

112: .seealso: [](ch:mfn), `MFNSolve()`
113: @*/
114: PetscErrorCode MFNSolveTranspose(MFN mfn,Vec b,Vec x)
115: {
116:   PetscFunctionBegin;
119:   PetscCheckSameComm(mfn,1,b,2);
121:   if (b!=x) PetscCheckSameComm(mfn,1,x,3);
122:   mfn->transpose_solve = PETSC_TRUE;
123:   if (!mfn->AT) PetscCall(MatCreateTranspose(mfn->A,&mfn->AT));
124:   PetscCall(MFNSolve_Private(mfn,b,x));
125:   PetscFunctionReturn(PETSC_SUCCESS);
126: }

128: /*@
129:    MFNGetIterationNumber - Gets the current iteration number. If the
130:    call to MFNSolve() is complete, then it returns the number of iterations
131:    carried out by the solution method.

133:    Not Collective

135:    Input Parameter:
136: .  mfn - the matrix function solver context

138:    Output Parameter:
139: .  its - number of iterations

141:    Note:
142:    During the i-th iteration this call returns i-1. If MFNSolve() is
143:    complete, then parameter "its" contains either the iteration number at
144:    which convergence was successfully reached, or failure was detected.
145:    Call MFNGetConvergedReason() to determine if the solver converged or
146:    failed and why.

148:    Level: intermediate

150: .seealso: [](ch:mfn), `MFNGetConvergedReason()`, `MFNSetTolerances()`
151: @*/
152: PetscErrorCode MFNGetIterationNumber(MFN mfn,PetscInt *its)
153: {
154:   PetscFunctionBegin;
156:   PetscAssertPointer(its,2);
157:   *its = mfn->its;
158:   PetscFunctionReturn(PETSC_SUCCESS);
159: }

161: /*@
162:    MFNGetConvergedReason - Gets the reason why the `MFNSolve()` iteration was
163:    stopped.

165:    Not Collective

167:    Input Parameter:
168: .  mfn - the matrix function solver context

170:    Output Parameter:
171: .  reason - negative value indicates diverged, positive value converged, see
172:    `MFNConvergedReason` for the possible values

174:    Options Database Key:
175: .  -mfn_converged_reason - print reason for convergence/divergence, and number of iterations

177:    Notes:
178:    If this routine is called before or doing the `MFNSolve()` the value of
179:    `MFN_CONVERGED_ITERATING` is returned.

181:    Basic solvers (e.g., unrestarted Krylov iterations) cannot determine if the
182:    computation is accurate up to the requested tolerance. In that case, the
183:    converged reason is set to `MFN_CONVERGED_ITS` if the requested number of steps
184:    (for instance, the `ncv` value in unrestarted Krylov methods) have been
185:    completed successfully.

187:    Level: intermediate

189: .seealso: [](ch:mfn), `MFNSetTolerances()`, `MFNSolve()`, `MFNConvergedReason`, `MFNSetErrorIfNotConverged()`
190: @*/
191: PetscErrorCode MFNGetConvergedReason(MFN mfn,MFNConvergedReason *reason)
192: {
193:   PetscFunctionBegin;
195:   PetscAssertPointer(reason,2);
196:   *reason = mfn->reason;
197:   PetscFunctionReturn(PETSC_SUCCESS);
198: }