Actual source code: lmeopts.c

slepc-main 2024-11-22
Report Typos and Errors
  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:    LME routines related to options that can be set via the command-line
 12:    or procedurally
 13: */

 15: #include <slepc/private/lmeimpl.h>
 16: #include <petscdraw.h>

 18: /*@C
 19:    LMEMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type
 20:    indicated by the user.

 22:    Collective

 24:    Input Parameters:
 25: +  lme      - the linear matrix equation context
 26: .  opt  - the command line option for this monitor
 27: .  name - the monitor type one is seeking
 28: -  ctx  - an optional user context for the monitor, or NULL

 30:    Level: developer

 32: .seealso: LMEMonitorSet()
 33: @*/
 34: PetscErrorCode LMEMonitorSetFromOptions(LME lme,const char opt[],const char name[],void *ctx)
 35: {
 36:   PetscErrorCode       (*mfunc)(LME,PetscInt,PetscReal,void*);
 37:   PetscErrorCode       (*cfunc)(PetscViewer,PetscViewerFormat,void*,PetscViewerAndFormat**);
 38:   PetscErrorCode       (*dfunc)(PetscViewerAndFormat**);
 39:   PetscViewerAndFormat *vf;
 40:   PetscViewer          viewer;
 41:   PetscViewerFormat    format;
 42:   PetscViewerType      vtype;
 43:   char                 key[PETSC_MAX_PATH_LEN];
 44:   PetscBool            flg;

 46:   PetscFunctionBegin;
 47:   PetscCall(PetscOptionsCreateViewer(PetscObjectComm((PetscObject)lme),((PetscObject)lme)->options,((PetscObject)lme)->prefix,opt,&viewer,&format,&flg));
 48:   if (!flg) PetscFunctionReturn(PETSC_SUCCESS);

 50:   PetscCall(PetscViewerGetType(viewer,&vtype));
 51:   PetscCall(SlepcMonitorMakeKey_Internal(name,vtype,format,key));
 52:   PetscCall(PetscFunctionListFind(LMEMonitorList,key,&mfunc));
 53:   PetscCheck(mfunc,PetscObjectComm((PetscObject)lme),PETSC_ERR_SUP,"Specified viewer and format not supported");
 54:   PetscCall(PetscFunctionListFind(LMEMonitorCreateList,key,&cfunc));
 55:   PetscCall(PetscFunctionListFind(LMEMonitorDestroyList,key,&dfunc));
 56:   if (!cfunc) cfunc = PetscViewerAndFormatCreate_Internal;
 57:   if (!dfunc) dfunc = PetscViewerAndFormatDestroy;

 59:   PetscCall((*cfunc)(viewer,format,ctx,&vf));
 60:   PetscCall(PetscViewerDestroy(&viewer));
 61:   PetscCall(LMEMonitorSet(lme,mfunc,vf,(PetscErrorCode(*)(void **))dfunc));
 62:   PetscFunctionReturn(PETSC_SUCCESS);
 63: }

 65: /*@
 66:    LMESetFromOptions - Sets LME options from the options database.
 67:    This routine must be called before LMESetUp() if the user is to be
 68:    allowed to set the solver type.

 70:    Collective

 72:    Input Parameters:
 73: .  lme - the linear matrix equation solver context

 75:    Notes:
 76:    To see all options, run your program with the -help option.

 78:    Level: beginner

 80: .seealso: LMESetOptionsPrefix()
 81: @*/
 82: PetscErrorCode LMESetFromOptions(LME lme)
 83: {
 84:   char           type[256];
 85:   PetscBool      set,flg,flg1,flg2;
 86:   PetscReal      r;
 87:   PetscInt       i;

 89:   PetscFunctionBegin;
 91:   PetscCall(LMERegisterAll());
 92:   PetscObjectOptionsBegin((PetscObject)lme);
 93:     PetscCall(PetscOptionsFList("-lme_type","Linear matrix equation","LMESetType",LMEList,(char*)(((PetscObject)lme)->type_name?((PetscObject)lme)->type_name:LMEKRYLOV),type,sizeof(type),&flg));
 94:     if (flg) PetscCall(LMESetType(lme,type));
 95:     else if (!((PetscObject)lme)->type_name) PetscCall(LMESetType(lme,LMEKRYLOV));

 97:     PetscCall(PetscOptionsBoolGroupBegin("-lme_lyapunov","Continuous-time Lyapunov equation","LMESetProblemType",&flg));
 98:     if (flg) PetscCall(LMESetProblemType(lme,LME_LYAPUNOV));
 99:     PetscCall(PetscOptionsBoolGroup("-lme_sylvester","Continuous-time Sylvester equation","LMESetProblemType",&flg));
100:     if (flg) PetscCall(LMESetProblemType(lme,LME_SYLVESTER));
101:     PetscCall(PetscOptionsBoolGroup("-lme_gen_lyapunov","Generalized Lyapunov equation","LMESetProblemType",&flg));
102:     if (flg) PetscCall(LMESetProblemType(lme,LME_GEN_LYAPUNOV));
103:     PetscCall(PetscOptionsBoolGroup("-lme_gen_sylvester","Generalized Sylvester equation","LMESetProblemType",&flg));
104:     if (flg) PetscCall(LMESetProblemType(lme,LME_GEN_SYLVESTER));
105:     PetscCall(PetscOptionsBoolGroup("-lme_dt_lyapunov","Discrete-time Lyapunov equation","LMESetProblemType",&flg));
106:     if (flg) PetscCall(LMESetProblemType(lme,LME_DT_LYAPUNOV));
107:     PetscCall(PetscOptionsBoolGroupEnd("-lme_stein","Stein equation","LMESetProblemType",&flg));
108:     if (flg) PetscCall(LMESetProblemType(lme,LME_STEIN));

110:     i = lme->max_it;
111:     PetscCall(PetscOptionsInt("-lme_max_it","Maximum number of iterations","LMESetTolerances",lme->max_it,&i,&flg1));
112:     if (!flg1) i = PETSC_DETERMINE;
113:     r = lme->tol;
114:     PetscCall(PetscOptionsReal("-lme_tol","Tolerance","LMESetTolerances",SlepcDefaultTol(lme->tol),&r,&flg2));
115:     if (flg1 || flg2) PetscCall(LMESetTolerances(lme,r,i));

117:     PetscCall(PetscOptionsInt("-lme_ncv","Number of basis vectors","LMESetDimensions",lme->ncv,&i,&flg));
118:     if (flg) PetscCall(LMESetDimensions(lme,i));

120:     PetscCall(PetscOptionsBool("-lme_error_if_not_converged","Generate error if solver does not converge","LMESetErrorIfNotConverged",lme->errorifnotconverged,&lme->errorifnotconverged,NULL));

122:     /* -----------------------------------------------------------------------*/
123:     /*
124:       Cancels all monitors hardwired into code before call to LMESetFromOptions()
125:     */
126:     PetscCall(PetscOptionsBool("-lme_monitor_cancel","Remove any hardwired monitor routines","LMEMonitorCancel",PETSC_FALSE,&flg,&set));
127:     if (set && flg) PetscCall(LMEMonitorCancel(lme));
128:     PetscCall(LMEMonitorSetFromOptions(lme,"-lme_monitor","error_estimate",NULL));

130:     /* -----------------------------------------------------------------------*/
131:     PetscCall(PetscOptionsName("-lme_view","Print detailed information on solver used","LMEView",&set));

133:     PetscTryTypeMethod(lme,setfromoptions,PetscOptionsObject);
134:     PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)lme,PetscOptionsObject));
135:   PetscOptionsEnd();

137:   if (!lme->V) PetscCall(LMEGetBV(lme,&lme->V));
138:   PetscCall(BVSetFromOptions(lme->V));
139:   PetscFunctionReturn(PETSC_SUCCESS);
140: }

142: /*@
143:    LMESetProblemType - Specifies the type of matrix equation to be solved.

145:    Logically Collective

147:    Input Parameters:
148: +  lme  - the linear matrix equation solver context
149: -  type - a known type of matrix equation

151:    Options Database Keys:
152: +  -lme_lyapunov - continuous-time Lyapunov equation A*X+X*A'=-C
153: .  -lme_sylvester - continuous-time Sylvester equation A*X+X*B=C
154: .  -lme_gen_lyapunov - generalized Lyapunov equation A*X*D'+D*X*A'=-C
155: .  -lme_gen_sylvester - generalized Sylvester equation A*X*E+D*X*B=C
156: .  -lme_dt_lyapunov - discrete-time Lyapunov equation A*X*A'-X=-C
157: -  -lme_stein - Stein equation A*X*E+X=C

159:    Notes:
160:    The coefficient matrices A, B, D, E must be provided via LMESetCoefficients(),
161:    but some of them are optional depending on the matrix equation.

163: .vb
164:                             equation              A    B    D    E
165:                           -----------------      ---  ---  ---  ---
166:        LME_LYAPUNOV        A*X+X*A'=-C           yes (A-t)  -    -
167:        LME_SYLVESTER       A*X+X*B=C             yes  yes   -    -
168:        LME_GEN_LYAPUNOV    A*X*D'+D*X*A'=-C      yes (A-t) yes (D-t)
169:        LME_GEN_SYLVESTER   A*X*E+D*X*B=C         yes  yes  yes  yes
170:        LME_DT_LYAPUNOV     A*X*A'-X=-C           yes   -    -  (A-t)
171:        LME_STEIN           A*X*E+X=C             yes   -    -   yes
172: .ve

174:    In the above table, the notation (A-t) means that this matrix need
175:    not be passed, but the user may choose to pass an explicit transpose
176:    of matrix A (for improved efficiency).

178:    Also note that some of the equation types impose restrictions on the
179:    properties of the coefficient matrices and possibly on the right-hand
180:    side C.

182:    Level: beginner

184: .seealso: LMESetCoefficients(), LMESetType(), LMEGetProblemType(), LMEProblemType
185: @*/
186: PetscErrorCode LMESetProblemType(LME lme,LMEProblemType type)
187: {
188:   PetscFunctionBegin;
191:   if (type == lme->problem_type) PetscFunctionReturn(PETSC_SUCCESS);
192:   switch (type) {
193:     case LME_LYAPUNOV:
194:     case LME_SYLVESTER:
195:     case LME_GEN_LYAPUNOV:
196:     case LME_GEN_SYLVESTER:
197:     case LME_DT_LYAPUNOV:
198:     case LME_STEIN:
199:       break;
200:     default:
201:       SETERRQ(PetscObjectComm((PetscObject)lme),PETSC_ERR_ARG_WRONG,"Unknown matrix equation type");
202:   }
203:   lme->problem_type = type;
204:   lme->setupcalled  = PETSC_FALSE;
205:   PetscFunctionReturn(PETSC_SUCCESS);
206: }

208: /*@
209:    LMEGetProblemType - Gets the matrix equation type from the LME object.

211:    Not Collective

213:    Input Parameter:
214: .  lme - the linear matrix equation solver context

216:    Output Parameter:
217: .  type - name of LME problem type

219:    Level: intermediate

221: .seealso: LMESetProblemType(), LMEProblemType
222: @*/
223: PetscErrorCode LMEGetProblemType(LME lme,LMEProblemType *type)
224: {
225:   PetscFunctionBegin;
227:   PetscAssertPointer(type,2);
228:   *type = lme->problem_type;
229:   PetscFunctionReturn(PETSC_SUCCESS);
230: }

232: /*@
233:    LMEGetTolerances - Gets the tolerance and maximum iteration count used
234:    by the LME convergence tests.

236:    Not Collective

238:    Input Parameter:
239: .  lme - the linear matrix equation solver context

241:    Output Parameters:
242: +  tol - the convergence tolerance
243: -  maxits - maximum number of iterations

245:    Notes:
246:    The user can specify NULL for any parameter that is not needed.

248:    Level: intermediate

250: .seealso: LMESetTolerances()
251: @*/
252: PetscErrorCode LMEGetTolerances(LME lme,PetscReal *tol,PetscInt *maxits)
253: {
254:   PetscFunctionBegin;
256:   if (tol)    *tol    = lme->tol;
257:   if (maxits) *maxits = lme->max_it;
258:   PetscFunctionReturn(PETSC_SUCCESS);
259: }

261: /*@
262:    LMESetTolerances - Sets the tolerance and maximum iteration count used
263:    by the LME convergence tests.

265:    Logically Collective

267:    Input Parameters:
268: +  lme - the linear matrix equation solver context
269: .  tol - the convergence tolerance
270: -  maxits - maximum number of iterations to use

272:    Options Database Keys:
273: +  -lme_tol <tol> - Sets the convergence tolerance
274: -  -lme_max_it <maxits> - Sets the maximum number of iterations allowed

276:    Notes:
277:    Use PETSC_CURRENT to retain the current value of any of the parameters.
278:    Use PETSC_DETERMINE for either argument to assign a default value computed
279:    internally (may be different in each solver).
280:    For maxits use PETSC_UMLIMITED to indicate there is no upper bound on this value.

282:    Level: intermediate

284: .seealso: LMEGetTolerances()
285: @*/
286: PetscErrorCode LMESetTolerances(LME lme,PetscReal tol,PetscInt maxits)
287: {
288:   PetscFunctionBegin;
292:   if (tol == (PetscReal)PETSC_DETERMINE) {
293:     lme->tol = PETSC_DETERMINE;
294:     lme->setupcalled = 0;
295:   } else if (tol != (PetscReal)PETSC_CURRENT) {
296:     PetscCheck(tol>0.0,PetscObjectComm((PetscObject)lme),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of tol. Must be > 0");
297:     lme->tol = tol;
298:   }
299:   if (maxits == PETSC_DETERMINE) {
300:     lme->max_it = PETSC_DETERMINE;
301:     lme->setupcalled = 0;
302:   } else if (maxits == PETSC_UNLIMITED) {
303:     lme->max_it = PETSC_INT_MAX;
304:   } else if (maxits != PETSC_CURRENT) {
305:     PetscCheck(maxits>0,PetscObjectComm((PetscObject)lme),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of maxits. Must be > 0");
306:     lme->max_it = maxits;
307:   }
308:   PetscFunctionReturn(PETSC_SUCCESS);
309: }

311: /*@
312:    LMEGetDimensions - Gets the dimension of the subspace used by the solver.

314:    Not Collective

316:    Input Parameter:
317: .  lme - the linear matrix equation solver context

319:    Output Parameter:
320: .  ncv - the maximum dimension of the subspace to be used by the solver

322:    Level: intermediate

324: .seealso: LMESetDimensions()
325: @*/
326: PetscErrorCode LMEGetDimensions(LME lme,PetscInt *ncv)
327: {
328:   PetscFunctionBegin;
330:   PetscAssertPointer(ncv,2);
331:   *ncv = lme->ncv;
332:   PetscFunctionReturn(PETSC_SUCCESS);
333: }

335: /*@
336:    LMESetDimensions - Sets the dimension of the subspace to be used by the solver.

338:    Logically Collective

340:    Input Parameters:
341: +  lme - the linear matrix equation solver context
342: -  ncv - the maximum dimension of the subspace to be used by the solver

344:    Options Database Keys:
345: .  -lme_ncv <ncv> - Sets the dimension of the subspace

347:    Notes:
348:    Use PETSC_DETERMINE for ncv to assign a reasonably good value, which is
349:    dependent on the solution method.

351:    Level: intermediate

353: .seealso: LMEGetDimensions()
354: @*/
355: PetscErrorCode LMESetDimensions(LME lme,PetscInt ncv)
356: {
357:   PetscFunctionBegin;
360:   if (ncv == PETSC_DECIDE || ncv == PETSC_DEFAULT) {
361:     lme->ncv = PETSC_DETERMINE;
362:   } else {
363:     PetscCheck(ncv>0,PetscObjectComm((PetscObject)lme),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of ncv. Must be > 0");
364:     lme->ncv = ncv;
365:   }
366:   lme->setupcalled = 0;
367:   PetscFunctionReturn(PETSC_SUCCESS);
368: }

370: /*@
371:    LMESetErrorIfNotConverged - Causes LMESolve() to generate an error if the
372:    solver has not converged.

374:    Logically Collective

376:    Input Parameters:
377: +  lme - the linear matrix equation solver context
378: -  flg - PETSC_TRUE indicates you want the error generated

380:    Options Database Keys:
381: .  -lme_error_if_not_converged - this takes an optional truth value (0/1/no/yes/true/false)

383:    Level: intermediate

385:    Note:
386:    Normally SLEPc continues if the solver fails to converge, you can call
387:    LMEGetConvergedReason() after a LMESolve() to determine if it has converged.

389: .seealso: LMEGetErrorIfNotConverged()
390: @*/
391: PetscErrorCode LMESetErrorIfNotConverged(LME lme,PetscBool flg)
392: {
393:   PetscFunctionBegin;
396:   lme->errorifnotconverged = flg;
397:   PetscFunctionReturn(PETSC_SUCCESS);
398: }

400: /*@
401:    LMEGetErrorIfNotConverged - Return a flag indicating whether LMESolve() will
402:    generate an error if the solver does not converge.

404:    Not Collective

406:    Input Parameter:
407: .  lme - the linear matrix equation solver context

409:    Output Parameter:
410: .  flag - PETSC_TRUE if it will generate an error, else PETSC_FALSE

412:    Level: intermediate

414: .seealso: LMESetErrorIfNotConverged()
415: @*/
416: PetscErrorCode LMEGetErrorIfNotConverged(LME lme,PetscBool *flag)
417: {
418:   PetscFunctionBegin;
420:   PetscAssertPointer(flag,2);
421:   *flag = lme->errorifnotconverged;
422:   PetscFunctionReturn(PETSC_SUCCESS);
423: }

425: /*@
426:    LMESetOptionsPrefix - Sets the prefix used for searching for all
427:    LME options in the database.

429:    Logically Collective

431:    Input Parameters:
432: +  lme - the linear matrix equation solver context
433: -  prefix - the prefix string to prepend to all LME option requests

435:    Notes:
436:    A hyphen (-) must NOT be given at the beginning of the prefix name.
437:    The first character of all runtime options is AUTOMATICALLY the
438:    hyphen.

440:    For example, to distinguish between the runtime options for two
441:    different LME contexts, one could call
442: .vb
443:       LMESetOptionsPrefix(lme1,"fun1_")
444:       LMESetOptionsPrefix(lme2,"fun2_")
445: .ve

447:    Level: advanced

449: .seealso: LMEAppendOptionsPrefix(), LMEGetOptionsPrefix()
450: @*/
451: PetscErrorCode LMESetOptionsPrefix(LME lme,const char *prefix)
452: {
453:   PetscFunctionBegin;
455:   if (!lme->V) PetscCall(LMEGetBV(lme,&lme->V));
456:   PetscCall(BVSetOptionsPrefix(lme->V,prefix));
457:   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)lme,prefix));
458:   PetscFunctionReturn(PETSC_SUCCESS);
459: }

461: /*@
462:    LMEAppendOptionsPrefix - Appends to the prefix used for searching for all
463:    LME options in the database.

465:    Logically Collective

467:    Input Parameters:
468: +  lme - the linear matrix equation solver context
469: -  prefix - the prefix string to prepend to all LME option requests

471:    Notes:
472:    A hyphen (-) must NOT be given at the beginning of the prefix name.
473:    The first character of all runtime options is AUTOMATICALLY the hyphen.

475:    Level: advanced

477: .seealso: LMESetOptionsPrefix(), LMEGetOptionsPrefix()
478: @*/
479: PetscErrorCode LMEAppendOptionsPrefix(LME lme,const char *prefix)
480: {
481:   PetscFunctionBegin;
483:   if (!lme->V) PetscCall(LMEGetBV(lme,&lme->V));
484:   PetscCall(BVAppendOptionsPrefix(lme->V,prefix));
485:   PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)lme,prefix));
486:   PetscFunctionReturn(PETSC_SUCCESS);
487: }

489: /*@
490:    LMEGetOptionsPrefix - Gets the prefix used for searching for all
491:    LME options in the database.

493:    Not Collective

495:    Input Parameters:
496: .  lme - the linear matrix equation solver context

498:    Output Parameters:
499: .  prefix - pointer to the prefix string used is returned

501:    Note:
502:    On the Fortran side, the user should pass in a string 'prefix' of
503:    sufficient length to hold the prefix.

505:    Level: advanced

507: .seealso: LMESetOptionsPrefix(), LMEAppendOptionsPrefix()
508: @*/
509: PetscErrorCode LMEGetOptionsPrefix(LME lme,const char *prefix[])
510: {
511:   PetscFunctionBegin;
513:   PetscAssertPointer(prefix,2);
514:   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)lme,prefix));
515:   PetscFunctionReturn(PETSC_SUCCESS);
516: }