Actual source code: lmeopts.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:    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 solver 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: [](ch:lme), `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 configure the solver.

 70:    Collective

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

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

 78:    Level: beginner

 80: .seealso: [](ch:lme), `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
153: .  -lme_sylvester     - Sylvester equation
154: .  -lme_gen_lyapunov  - generalized Lyapunov equation
155: .  -lme_gen_sylvester - generalized Sylvester equation
156: .  -lme_dt_lyapunov   - discrete-time Lyapunov equation
157: -  -lme_stein         - Stein equation

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:    Problem Type             | Equation         |`LMEProblemType`   | $A$ | $B$ | $D$ | $E$
164:    -------------------------|------------------|-------------------|-----|-----|-----|-----
165:    Continuous-Time Lyapunov | $AX+XA^*=-C$     |`LME_LYAPUNOV`     | yes |$A^*$|  -  |  -
166:    Sylvester                | $AX+XB=C$        |`LME_SYLVESTER`    | yes | yes |  -  |  -
167:    Generalized Lyapunov     | $AXD^*+DXA^*=-C$ |`LME_GEN_LYAPUNOV` | yes |$A^*$| yes |$D^*$
168:    Generalized Sylvester    | $AXE+DXB=C$      |`LME_GEN_SYLVESTER`| yes | yes | yes | yes
169:    Discrete-Time Lyapunov   | $AXA^*-X=-C$     |`LME_DT_LYAPUNOV`  | yes |  -  |  -  |$A^*$
170:    Stein                    | $AXE-X=-C$       |`LME_STEIN`        | yes |  -  |  -  | yes

172:    In the above table, the notation $A^*$ means that this matrix need
173:    not be passed, but the user may choose to pass an explicit transpose
174:    of matrix $A$ (for improved efficiency).

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

180:    Level: beginner

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

206: /*@
207:    LMEGetProblemType - Gets the matrix equation type from the `LME` object.

209:    Not Collective

211:    Input Parameter:
212: .  lme - the linear matrix equation solver context

214:    Output Parameter:
215: .  type - the problem type

217:    Level: beginner

219: .seealso: [](ch:lme), `LMESetProblemType()`, `LMEProblemType`
220: @*/
221: PetscErrorCode LMEGetProblemType(LME lme,LMEProblemType *type)
222: {
223:   PetscFunctionBegin;
225:   PetscAssertPointer(type,2);
226:   *type = lme->problem_type;
227:   PetscFunctionReturn(PETSC_SUCCESS);
228: }

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

234:    Not Collective

236:    Input Parameter:
237: .  lme - the linear matrix equation solver context

239:    Output Parameters:
240: +  tol - the convergence tolerance
241: -  maxits - maximum number of iterations

243:    Note:
244:    The user can specify `NULL` for any parameter that is not needed.

246:    Level: intermediate

248: .seealso: [](ch:lme), `LMESetTolerances()`
249: @*/
250: PetscErrorCode LMEGetTolerances(LME lme,PetscReal *tol,PetscInt *maxits)
251: {
252:   PetscFunctionBegin;
254:   if (tol)    *tol    = lme->tol;
255:   if (maxits) *maxits = lme->max_it;
256:   PetscFunctionReturn(PETSC_SUCCESS);
257: }

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

263:    Logically Collective

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

270:    Options Database Keys:
271: +  -lme_tol \<tol\>       - sets the convergence tolerance
272: -  -lme_max_it \<maxits\> - sets the maximum number of iterations allowed

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

280:    Level: intermediate

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

309: /*@
310:    LMEGetDimensions - Gets the dimension of the subspace used by the solver.

312:    Not Collective

314:    Input Parameter:
315: .  lme - the linear matrix equation solver context

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

320:    Level: intermediate

322: .seealso: [](ch:lme), `LMESetDimensions()`
323: @*/
324: PetscErrorCode LMEGetDimensions(LME lme,PetscInt *ncv)
325: {
326:   PetscFunctionBegin;
328:   PetscAssertPointer(ncv,2);
329:   *ncv = lme->ncv;
330:   PetscFunctionReturn(PETSC_SUCCESS);
331: }

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

336:    Logically Collective

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

342:    Options Database Key:
343: .  -lme_ncv \<ncv\> - sets the dimension of the subspace

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

349:    Level: intermediate

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

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

372:    Logically Collective

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

378:    Options Database Key:
379: .  -lme_error_if_not_converged - generate an error and stop the program

381:    Note:
382:    Normally SLEPc continues if the solver fails to converge, you can call
383:    `LMEGetConvergedReason()` after a `LMESolve()` to determine if it has converged.

385:    Level: intermediate

387: .seealso: [](ch:lme), `LMEGetErrorIfNotConverged()`
388: @*/
389: PetscErrorCode LMESetErrorIfNotConverged(LME lme,PetscBool flg)
390: {
391:   PetscFunctionBegin;
394:   lme->errorifnotconverged = flg;
395:   PetscFunctionReturn(PETSC_SUCCESS);
396: }

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

402:    Not Collective

404:    Input Parameter:
405: .  lme - the linear matrix equation solver context

407:    Output Parameter:
408: .  flag - `PETSC_TRUE` if it will generate an error, else `PETSC_FALSE`

410:    Level: intermediate

412: .seealso: [](ch:lme), `LMESetErrorIfNotConverged()`
413: @*/
414: PetscErrorCode LMEGetErrorIfNotConverged(LME lme,PetscBool *flag)
415: {
416:   PetscFunctionBegin;
418:   PetscAssertPointer(flag,2);
419:   *flag = lme->errorifnotconverged;
420:   PetscFunctionReturn(PETSC_SUCCESS);
421: }

423: /*@
424:    LMESetOptionsPrefix - Sets the prefix used for searching for all
425:    `LME` options in the database.

427:    Logically Collective

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

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

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

445:    Level: advanced

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

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

463:    Logically Collective

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

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

473:    Level: advanced

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

487: /*@
488:    LMEGetOptionsPrefix - Gets the prefix used for searching for all
489:    `LME` options in the database.

491:    Not Collective

493:    Input Parameter:
494: .  lme - the linear matrix equation solver context

496:    Output Parameter:
497: .  prefix - pointer to the prefix string used is returned

499:    Level: advanced

501: .seealso: [](ch:lme), `LMESetOptionsPrefix()`, `LMEAppendOptionsPrefix()`
502: @*/
503: PetscErrorCode LMEGetOptionsPrefix(LME lme,const char *prefix[])
504: {
505:   PetscFunctionBegin;
507:   PetscAssertPointer(prefix,2);
508:   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)lme,prefix));
509:   PetscFunctionReturn(PETSC_SUCCESS);
510: }