| Line | Branch | Exec | Source |
|---|---|---|---|
| 1 | /* | ||
| 2 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 3 | SLEPc - Scalable Library for Eigenvalue Problem Computations | ||
| 4 | Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain | ||
| 5 | |||
| 6 | This file is part of SLEPc. | ||
| 7 | SLEPc is distributed under a 2-clause BSD license (see LICENSE). | ||
| 8 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 9 | */ | ||
| 10 | /* | ||
| 11 | NEP routines related to options that can be set via the command-line | ||
| 12 | or procedurally | ||
| 13 | */ | ||
| 14 | |||
| 15 | #include <slepc/private/nepimpl.h> /*I "slepcnep.h" I*/ | ||
| 16 | #include <petscdraw.h> | ||
| 17 | |||
| 18 | /*@C | ||
| 19 | NEPMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type | ||
| 20 | indicated by the user. | ||
| 21 | |||
| 22 | Collective | ||
| 23 | |||
| 24 | Input Parameters: | ||
| 25 | + nep - the nonlinear eigensolver 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` | ||
| 29 | - trackall - whether this monitor tracks all eigenvalues or not | ||
| 30 | |||
| 31 | Level: developer | ||
| 32 | |||
| 33 | .seealso: [](ch:nep), `NEPMonitorSet()`, `NEPSetTrackAll()` | ||
| 34 | @*/ | ||
| 35 | 3927 | PetscErrorCode NEPMonitorSetFromOptions(NEP nep,const char opt[],const char name[],void *ctx,PetscBool trackall) | |
| 36 | { | ||
| 37 | 3927 | PetscErrorCode (*mfunc)(NEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*); | |
| 38 | 3927 | PetscErrorCode (*cfunc)(PetscViewer,PetscViewerFormat,void*,PetscViewerAndFormat**); | |
| 39 | 3927 | PetscErrorCode (*dfunc)(PetscViewerAndFormat**); | |
| 40 | 3927 | PetscViewerAndFormat *vf; | |
| 41 | 3927 | PetscViewer viewer; | |
| 42 | 3927 | PetscViewerFormat format; | |
| 43 | 3927 | PetscViewerType vtype; | |
| 44 | 3927 | char key[PETSC_MAX_PATH_LEN]; | |
| 45 | 3927 | PetscBool flg; | |
| 46 | |||
| 47 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
3927 | PetscFunctionBegin; |
| 48 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
3927 | PetscCall(PetscOptionsCreateViewer(PetscObjectComm((PetscObject)nep),((PetscObject)nep)->options,((PetscObject)nep)->prefix,opt,&viewer,&format,&flg)); |
| 49 |
8/14✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✓ Branch 13 taken 2 times.
|
3927 | if (!flg) PetscFunctionReturn(PETSC_SUCCESS); |
| 50 | |||
| 51 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
54 | PetscCall(PetscViewerGetType(viewer,&vtype)); |
| 52 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
54 | PetscCall(SlepcMonitorMakeKey_Internal(name,vtype,format,key)); |
| 53 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
54 | PetscCall(PetscFunctionListFind(NEPMonitorList,key,&mfunc)); |
| 54 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
54 | PetscCheck(mfunc,PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"Specified viewer and format not supported"); |
| 55 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
54 | PetscCall(PetscFunctionListFind(NEPMonitorCreateList,key,&cfunc)); |
| 56 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
54 | PetscCall(PetscFunctionListFind(NEPMonitorDestroyList,key,&dfunc)); |
| 57 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
54 | if (!cfunc) cfunc = PetscViewerAndFormatCreate_Internal; |
| 58 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
54 | if (!dfunc) dfunc = PetscViewerAndFormatDestroy; |
| 59 | |||
| 60 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
54 | PetscCall((*cfunc)(viewer,format,ctx,&vf)); |
| 61 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
54 | PetscCall(PetscViewerDestroy(&viewer)); |
| 62 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
54 | PetscCall(NEPMonitorSet(nep,mfunc,vf,(PetscCtxDestroyFn*)dfunc)); |
| 63 |
6/8✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
54 | if (trackall) PetscCall(NEPSetTrackAll(nep,PETSC_TRUE)); |
| 64 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
9 | PetscFunctionReturn(PETSC_SUCCESS); |
| 65 | } | ||
| 66 | |||
| 67 | /*@ | ||
| 68 | NEPSetFromOptions - Sets `NEP` options from the options database. | ||
| 69 | This routine must be called before `NEPSetUp()` if the user is to be | ||
| 70 | allowed to configure the solver. | ||
| 71 | |||
| 72 | Collective | ||
| 73 | |||
| 74 | Input Parameter: | ||
| 75 | . nep - the nonlinear eigensolver context | ||
| 76 | |||
| 77 | Note: | ||
| 78 | To see all options, run your program with the `-help` option. | ||
| 79 | |||
| 80 | Level: beginner | ||
| 81 | |||
| 82 | .seealso: [](ch:nep), `NEPSetOptionsPrefix()` | ||
| 83 | @*/ | ||
| 84 | 1309 | PetscErrorCode NEPSetFromOptions(NEP nep) | |
| 85 | { | ||
| 86 | 1309 | char type[256]; | |
| 87 | 1309 | PetscBool set,flg,flg1,flg2,flg3,flg4,flg5,bval; | |
| 88 | 1309 | PetscReal r; | |
| 89 | 1309 | PetscScalar s; | |
| 90 | 1309 | PetscInt i,j,k; | |
| 91 | 1309 | NEPRefine refine; | |
| 92 | 1309 | NEPRefineScheme scheme; | |
| 93 | |||
| 94 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
1309 | PetscFunctionBegin; |
| 95 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
1309 | PetscValidHeaderSpecific(nep,NEP_CLASSID,1); |
| 96 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1309 | PetscCall(NEPRegisterAll()); |
| 97 |
8/10✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 10 times.
✓ Branch 5 taken 8 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
✓ Branch 8 taken 2 times.
✓ Branch 9 taken 2 times.
|
3927 | PetscObjectOptionsBegin((PetscObject)nep); |
| 98 |
7/8✓ Branch 0 taken 10 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
1759 | PetscCall(PetscOptionsFList("-nep_type","Nonlinear eigensolver method","NEPSetType",NEPList,(char*)(((PetscObject)nep)->type_name?((PetscObject)nep)->type_name:NEPRII),type,sizeof(type),&flg)); |
| 99 |
6/8✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
1309 | if (flg) PetscCall(NEPSetType(nep,type)); |
| 100 |
6/8✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
369 | else if (!((PetscObject)nep)->type_name) PetscCall(NEPSetType(nep,NEPRII)); |
| 101 | |||
| 102 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1309 | PetscCall(PetscOptionsBoolGroupBegin("-nep_general","General nonlinear eigenvalue problem","NEPSetProblemType",&flg)); |
| 103 |
6/8✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
1309 | if (flg) PetscCall(NEPSetProblemType(nep,NEP_GENERAL)); |
| 104 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1309 | PetscCall(PetscOptionsBoolGroupEnd("-nep_rational","Rational eigenvalue problem","NEPSetProblemType",&flg)); |
| 105 |
1/8✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
1309 | if (flg) PetscCall(NEPSetProblemType(nep,NEP_RATIONAL)); |
| 106 | |||
| 107 | 1309 | refine = nep->refine; | |
| 108 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1309 | PetscCall(PetscOptionsEnum("-nep_refine","Iterative refinement method","NEPSetRefine",NEPRefineTypes,(PetscEnum)refine,(PetscEnum*)&refine,&flg1)); |
| 109 | 1309 | i = nep->npart; | |
| 110 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1309 | PetscCall(PetscOptionsInt("-nep_refine_partitions","Number of partitions of the communicator for iterative refinement","NEPSetRefine",nep->npart,&i,&flg2)); |
| 111 | 1309 | r = nep->rtol; | |
| 112 |
7/8✓ Branch 0 taken 10 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
1325 | PetscCall(PetscOptionsReal("-nep_refine_tol","Tolerance for iterative refinement","NEPSetRefine",nep->rtol==(PetscReal)PETSC_DETERMINE?SLEPC_DEFAULT_TOL/1000:nep->rtol,&r,&flg3)); |
| 113 | 1309 | j = nep->rits; | |
| 114 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1309 | PetscCall(PetscOptionsInt("-nep_refine_its","Maximum number of iterations for iterative refinement","NEPSetRefine",nep->rits,&j,&flg4)); |
| 115 | 1309 | scheme = nep->scheme; | |
| 116 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1309 | PetscCall(PetscOptionsEnum("-nep_refine_scheme","Scheme used for linear systems within iterative refinement","NEPSetRefine",NEPRefineSchemes,(PetscEnum)scheme,(PetscEnum*)&scheme,&flg5)); |
| 117 |
10/16✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 10 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 10 times.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 10 times.
✓ Branch 10 taken 2 times.
✓ Branch 11 taken 8 times.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
|
1309 | if (flg1 || flg2 || flg3 || flg4 || flg5) PetscCall(NEPSetRefine(nep,refine,i,r,j,scheme)); |
| 118 | |||
| 119 | 1309 | i = nep->max_it; | |
| 120 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1309 | PetscCall(PetscOptionsInt("-nep_max_it","Maximum number of iterations","NEPSetTolerances",nep->max_it,&i,&flg1)); |
| 121 | 1309 | r = nep->tol; | |
| 122 |
7/8✓ Branch 0 taken 10 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
2102 | PetscCall(PetscOptionsReal("-nep_tol","Tolerance","NEPSetTolerances",SlepcDefaultTol(nep->tol),&r,&flg2)); |
| 123 |
7/10✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 10 times.
✓ Branch 4 taken 2 times.
✓ Branch 5 taken 8 times.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
|
1309 | if (flg1 || flg2) PetscCall(NEPSetTolerances(nep,r,i)); |
| 124 | |||
| 125 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1309 | PetscCall(PetscOptionsBoolGroupBegin("-nep_conv_rel","Relative error convergence test","NEPSetConvergenceTest",&flg)); |
| 126 |
1/8✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
1309 | if (flg) PetscCall(NEPSetConvergenceTest(nep,NEP_CONV_REL)); |
| 127 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1309 | PetscCall(PetscOptionsBoolGroup("-nep_conv_norm","Convergence test relative to the matrix norms","NEPSetConvergenceTest",&flg)); |
| 128 |
6/8✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
1309 | if (flg) PetscCall(NEPSetConvergenceTest(nep,NEP_CONV_NORM)); |
| 129 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1309 | PetscCall(PetscOptionsBoolGroup("-nep_conv_abs","Absolute error convergence test","NEPSetConvergenceTest",&flg)); |
| 130 |
1/8✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
1309 | if (flg) PetscCall(NEPSetConvergenceTest(nep,NEP_CONV_ABS)); |
| 131 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1309 | PetscCall(PetscOptionsBoolGroupEnd("-nep_conv_user","User-defined convergence test","NEPSetConvergenceTest",&flg)); |
| 132 |
1/8✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
1309 | if (flg) PetscCall(NEPSetConvergenceTest(nep,NEP_CONV_USER)); |
| 133 | |||
| 134 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1309 | PetscCall(PetscOptionsBoolGroupBegin("-nep_stop_basic","Stop iteration if all eigenvalues converged or max_it reached","NEPSetStoppingTest",&flg)); |
| 135 |
1/8✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
1309 | if (flg) PetscCall(NEPSetStoppingTest(nep,NEP_STOP_BASIC)); |
| 136 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1309 | PetscCall(PetscOptionsBoolGroupEnd("-nep_stop_user","User-defined stopping test","NEPSetStoppingTest",&flg)); |
| 137 |
1/8✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
1309 | if (flg) PetscCall(NEPSetStoppingTest(nep,NEP_STOP_USER)); |
| 138 | |||
| 139 | 1309 | i = nep->nev; | |
| 140 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1309 | PetscCall(PetscOptionsInt("-nep_nev","Number of eigenvalues to compute","NEPSetDimensions",nep->nev,&i,&flg1)); |
| 141 | 1309 | j = nep->ncv; | |
| 142 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1309 | PetscCall(PetscOptionsInt("-nep_ncv","Number of basis vectors","NEPSetDimensions",nep->ncv,&j,&flg2)); |
| 143 | 1309 | k = nep->mpd; | |
| 144 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1309 | PetscCall(PetscOptionsInt("-nep_mpd","Maximum dimension of projected problem","NEPSetDimensions",nep->mpd,&k,&flg3)); |
| 145 |
9/12✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 10 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 10 times.
✓ Branch 6 taken 2 times.
✓ Branch 7 taken 8 times.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
1309 | if (flg1 || flg2 || flg3) PetscCall(NEPSetDimensions(nep,i,j,k)); |
| 146 | |||
| 147 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1309 | PetscCall(PetscOptionsBoolGroupBegin("-nep_largest_magnitude","Compute largest eigenvalues in magnitude","NEPSetWhichEigenpairs",&flg)); |
| 148 |
1/8✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
1309 | if (flg) PetscCall(NEPSetWhichEigenpairs(nep,NEP_LARGEST_MAGNITUDE)); |
| 149 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1309 | PetscCall(PetscOptionsBoolGroup("-nep_smallest_magnitude","Compute smallest eigenvalues in magnitude","NEPSetWhichEigenpairs",&flg)); |
| 150 |
1/8✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
1309 | if (flg) PetscCall(NEPSetWhichEigenpairs(nep,NEP_SMALLEST_MAGNITUDE)); |
| 151 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1309 | PetscCall(PetscOptionsBoolGroup("-nep_largest_real","Compute eigenvalues with largest real parts","NEPSetWhichEigenpairs",&flg)); |
| 152 |
1/8✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
1309 | if (flg) PetscCall(NEPSetWhichEigenpairs(nep,NEP_LARGEST_REAL)); |
| 153 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1309 | PetscCall(PetscOptionsBoolGroup("-nep_smallest_real","Compute eigenvalues with smallest real parts","NEPSetWhichEigenpairs",&flg)); |
| 154 |
1/8✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
1309 | if (flg) PetscCall(NEPSetWhichEigenpairs(nep,NEP_SMALLEST_REAL)); |
| 155 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1309 | PetscCall(PetscOptionsBoolGroup("-nep_largest_imaginary","Compute eigenvalues with largest imaginary parts","NEPSetWhichEigenpairs",&flg)); |
| 156 |
1/8✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
1309 | if (flg) PetscCall(NEPSetWhichEigenpairs(nep,NEP_LARGEST_IMAGINARY)); |
| 157 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1309 | PetscCall(PetscOptionsBoolGroup("-nep_smallest_imaginary","Compute eigenvalues with smallest imaginary parts","NEPSetWhichEigenpairs",&flg)); |
| 158 |
1/8✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
1309 | if (flg) PetscCall(NEPSetWhichEigenpairs(nep,NEP_SMALLEST_IMAGINARY)); |
| 159 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1309 | PetscCall(PetscOptionsBoolGroup("-nep_target_magnitude","Compute eigenvalues closest to target","NEPSetWhichEigenpairs",&flg)); |
| 160 |
1/8✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
1309 | if (flg) PetscCall(NEPSetWhichEigenpairs(nep,NEP_TARGET_MAGNITUDE)); |
| 161 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1309 | PetscCall(PetscOptionsBoolGroup("-nep_target_real","Compute eigenvalues with real parts closest to target","NEPSetWhichEigenpairs",&flg)); |
| 162 |
6/8✓ Branch 0 taken 5 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 1 times.
✓ Branch 3 taken 4 times.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 1 times.
|
1309 | if (flg) PetscCall(NEPSetWhichEigenpairs(nep,NEP_TARGET_REAL)); |
| 163 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1309 | PetscCall(PetscOptionsBoolGroup("-nep_target_imaginary","Compute eigenvalues with imaginary parts closest to target","NEPSetWhichEigenpairs",&flg)); |
| 164 |
1/8✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
1309 | if (flg) PetscCall(NEPSetWhichEigenpairs(nep,NEP_TARGET_IMAGINARY)); |
| 165 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1309 | PetscCall(PetscOptionsBoolGroup("-nep_all","Compute all eigenvalues in a region","NEPSetWhichEigenpairs",&flg)); |
| 166 |
1/8✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
1309 | if (flg) PetscCall(NEPSetWhichEigenpairs(nep,NEP_ALL)); |
| 167 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1309 | PetscCall(PetscOptionsBoolGroupEnd("-nep_which_user","Select the user-defined selection criterion","NEPSetWhichEigenpairs",&flg)); |
| 168 |
1/8✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
1309 | if (flg) PetscCall(NEPSetWhichEigenpairs(nep,NEP_WHICH_USER)); |
| 169 | |||
| 170 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1309 | PetscCall(PetscOptionsScalar("-nep_target","Value of the target","NEPSetTarget",nep->target,&s,&flg)); |
| 171 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
1309 | if (flg) { |
| 172 |
6/8✓ Branch 0 taken 10 times.
✓ Branch 1 taken 5 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
727 | if (nep->which!=NEP_TARGET_REAL && nep->which!=NEP_TARGET_IMAGINARY) PetscCall(NEPSetWhichEigenpairs(nep,NEP_TARGET_MAGNITUDE)); |
| 173 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
727 | PetscCall(NEPSetTarget(nep,s)); |
| 174 | } | ||
| 175 | |||
| 176 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1309 | PetscCall(PetscOptionsBool("-nep_two_sided","Use two-sided variant (to compute left eigenvectors)","NEPSetTwoSided",nep->twosided,&bval,&flg)); |
| 177 |
6/8✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
1309 | if (flg) PetscCall(NEPSetTwoSided(nep,bval)); |
| 178 | |||
| 179 | /* -----------------------------------------------------------------------*/ | ||
| 180 | /* | ||
| 181 | Cancels all monitors hardwired into code before call to NEPSetFromOptions() | ||
| 182 | */ | ||
| 183 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1309 | PetscCall(PetscOptionsBool("-nep_monitor_cancel","Remove any hardwired monitor routines","NEPMonitorCancel",PETSC_FALSE,&flg,&set)); |
| 184 |
1/10✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
|
1309 | if (set && flg) PetscCall(NEPMonitorCancel(nep)); |
| 185 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1309 | PetscCall(NEPMonitorSetFromOptions(nep,"-nep_monitor","first_approximation",NULL,PETSC_FALSE)); |
| 186 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1309 | PetscCall(NEPMonitorSetFromOptions(nep,"-nep_monitor_all","all_approximations",NULL,PETSC_TRUE)); |
| 187 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1309 | PetscCall(NEPMonitorSetFromOptions(nep,"-nep_monitor_conv","convergence_history",NULL,PETSC_FALSE)); |
| 188 | |||
| 189 | /* -----------------------------------------------------------------------*/ | ||
| 190 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1309 | PetscCall(PetscOptionsName("-nep_view","Print detailed information on solver used","NEPView",&set)); |
| 191 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1309 | PetscCall(PetscOptionsName("-nep_view_vectors","View computed eigenvectors","NEPVectorsView",&set)); |
| 192 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1309 | PetscCall(PetscOptionsName("-nep_view_values","View computed eigenvalues","NEPValuesView",&set)); |
| 193 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1309 | PetscCall(PetscOptionsName("-nep_converged_reason","Print reason for convergence, and number of iterations","NEPConvergedReasonView",&set)); |
| 194 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1309 | PetscCall(PetscOptionsName("-nep_error_absolute","Print absolute errors of each eigenpair","NEPErrorView",&set)); |
| 195 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1309 | PetscCall(PetscOptionsName("-nep_error_relative","Print relative errors of each eigenpair","NEPErrorView",&set)); |
| 196 | |||
| 197 |
5/8✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
1309 | PetscTryTypeMethod(nep,setfromoptions,PetscOptionsObject); |
| 198 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1309 | PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)nep,PetscOptionsObject)); |
| 199 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1309 | PetscOptionsEnd(); |
| 200 | |||
| 201 |
6/8✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
1309 | if (!nep->V) PetscCall(NEPGetBV(nep,&nep->V)); |
| 202 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1309 | PetscCall(BVSetFromOptions(nep->V)); |
| 203 |
6/8✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
1309 | if (!nep->rg) PetscCall(NEPGetRG(nep,&nep->rg)); |
| 204 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1309 | PetscCall(RGSetFromOptions(nep->rg)); |
| 205 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
1309 | if (nep->useds) { |
| 206 |
6/8✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
1059 | if (!nep->ds) PetscCall(NEPGetDS(nep,&nep->ds)); |
| 207 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1059 | PetscCall(NEPSetDSType(nep)); |
| 208 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1059 | PetscCall(DSSetFromOptions(nep->ds)); |
| 209 | } | ||
| 210 |
5/8✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
1309 | if (!nep->refineksp) PetscCall(NEPRefineGetKSP(nep,&nep->refineksp)); |
| 211 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1309 | PetscCall(KSPSetFromOptions(nep->refineksp)); |
| 212 |
9/10✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 10 times.
✓ Branch 5 taken 8 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
✓ Branch 8 taken 2 times.
✓ Branch 9 taken 2 times.
|
3888 | if (nep->fui==NEP_USER_INTERFACE_SPLIT) for (i=0;i<nep->nt;i++) PetscCall(FNSetFromOptions(nep->f[i])); |
| 213 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
245 | PetscFunctionReturn(PETSC_SUCCESS); |
| 214 | } | ||
| 215 | |||
| 216 | /*@ | ||
| 217 | NEPGetTolerances - Gets the tolerance and maximum iteration count used | ||
| 218 | by the `NEP` convergence tests. | ||
| 219 | |||
| 220 | Not Collective | ||
| 221 | |||
| 222 | Input Parameter: | ||
| 223 | . nep - the nonlinear eigensolver context | ||
| 224 | |||
| 225 | Output Parameters: | ||
| 226 | + tol - the convergence tolerance | ||
| 227 | - maxits - maximum number of iterations | ||
| 228 | |||
| 229 | Notes: | ||
| 230 | The user can specify `NULL` for any parameter that is not needed. | ||
| 231 | |||
| 232 | Level: intermediate | ||
| 233 | |||
| 234 | .seealso: [](ch:nep), `NEPSetTolerances()` | ||
| 235 | @*/ | ||
| 236 | 87 | PetscErrorCode NEPGetTolerances(NEP nep,PetscReal *tol,PetscInt *maxits) | |
| 237 | { | ||
| 238 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
87 | PetscFunctionBegin; |
| 239 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
87 | PetscValidHeaderSpecific(nep,NEP_CLASSID,1); |
| 240 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
87 | if (tol) *tol = nep->tol; |
| 241 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
87 | if (maxits) *maxits = nep->max_it; |
| 242 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
87 | PetscFunctionReturn(PETSC_SUCCESS); |
| 243 | } | ||
| 244 | |||
| 245 | /*@ | ||
| 246 | NEPSetTolerances - Sets the tolerance and maximum iteration count used | ||
| 247 | by the `NEP` convergence tests. | ||
| 248 | |||
| 249 | Logically Collective | ||
| 250 | |||
| 251 | Input Parameters: | ||
| 252 | + nep - the nonlinear eigensolver context | ||
| 253 | . tol - the convergence tolerance | ||
| 254 | - maxits - maximum number of iterations to use | ||
| 255 | |||
| 256 | Options Database Keys: | ||
| 257 | + -nep_tol \<tol\> - sets the convergence tolerance | ||
| 258 | - -nep_max_it \<maxits\> - sets the maximum number of iterations allowed | ||
| 259 | |||
| 260 | Note: | ||
| 261 | Use `PETSC_CURRENT` to retain the current value of any of the parameters. | ||
| 262 | Use `PETSC_DETERMINE` for either argument to assign a default value computed | ||
| 263 | internally (may be different in each solver). | ||
| 264 | For `maxits` use `PETSC_UNLIMITED` to indicate there is no upper bound on this value. | ||
| 265 | |||
| 266 | Level: intermediate | ||
| 267 | |||
| 268 | .seealso: [](ch:nep), `NEPGetTolerances()` | ||
| 269 | @*/ | ||
| 270 | 686 | PetscErrorCode NEPSetTolerances(NEP nep,PetscReal tol,PetscInt maxits) | |
| 271 | { | ||
| 272 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
686 | PetscFunctionBegin; |
| 273 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
686 | PetscValidHeaderSpecific(nep,NEP_CLASSID,1); |
| 274 |
29/66✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
✗ Branch 16 not taken.
✓ Branch 17 taken 2 times.
✓ Branch 18 taken 2 times.
✗ Branch 19 not taken.
✗ Branch 20 not taken.
✓ Branch 21 taken 2 times.
✓ Branch 22 taken 2 times.
✗ Branch 23 not taken.
✗ Branch 24 not taken.
✓ Branch 25 taken 2 times.
✓ Branch 26 taken 2 times.
✗ Branch 27 not taken.
✗ Branch 28 not taken.
✓ Branch 29 taken 2 times.
✓ Branch 30 taken 2 times.
✗ Branch 31 not taken.
✗ Branch 32 not taken.
✓ Branch 33 taken 2 times.
✓ Branch 34 taken 2 times.
✗ Branch 35 not taken.
✓ Branch 36 taken 2 times.
✗ Branch 37 not taken.
✗ Branch 38 not taken.
✓ Branch 39 taken 2 times.
✗ Branch 40 not taken.
✗ Branch 41 not taken.
✗ Branch 42 not taken.
✓ Branch 43 taken 2 times.
✗ Branch 44 not taken.
✗ Branch 45 not taken.
✗ Branch 46 not taken.
✓ Branch 47 taken 2 times.
✓ Branch 48 taken 2 times.
✗ Branch 49 not taken.
✗ Branch 50 not taken.
✓ Branch 51 taken 2 times.
✓ Branch 52 taken 2 times.
✗ Branch 53 not taken.
✓ Branch 54 taken 2 times.
✗ Branch 55 not taken.
✗ Branch 56 not taken.
✓ Branch 57 taken 2 times.
✗ Branch 58 not taken.
✗ Branch 59 not taken.
✓ Branch 60 taken 2 times.
✗ Branch 61 not taken.
✗ Branch 62 not taken.
✓ Branch 63 taken 2 times.
✗ Branch 64 not taken.
✗ Branch 65 not taken.
|
686 | PetscValidLogicalCollectiveReal(nep,tol,2); |
| 275 |
27/62✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✓ Branch 13 taken 2 times.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
✓ Branch 16 taken 2 times.
✗ Branch 17 not taken.
✗ Branch 18 not taken.
✓ Branch 19 taken 2 times.
✓ Branch 20 taken 2 times.
✗ Branch 21 not taken.
✗ Branch 22 not taken.
✓ Branch 23 taken 2 times.
✓ Branch 24 taken 2 times.
✗ Branch 25 not taken.
✗ Branch 26 not taken.
✓ Branch 27 taken 2 times.
✓ Branch 28 taken 2 times.
✗ Branch 29 not taken.
✗ Branch 30 not taken.
✓ Branch 31 taken 2 times.
✓ Branch 32 taken 2 times.
✗ Branch 33 not taken.
✓ Branch 34 taken 2 times.
✗ Branch 35 not taken.
✗ Branch 36 not taken.
✓ Branch 37 taken 2 times.
✗ Branch 38 not taken.
✗ Branch 39 not taken.
✗ Branch 40 not taken.
✓ Branch 41 taken 2 times.
✗ Branch 42 not taken.
✗ Branch 43 not taken.
✗ Branch 44 not taken.
✓ Branch 45 taken 2 times.
✓ Branch 46 taken 2 times.
✗ Branch 47 not taken.
✗ Branch 48 not taken.
✓ Branch 49 taken 2 times.
✓ Branch 50 taken 2 times.
✗ Branch 51 not taken.
✓ Branch 52 taken 2 times.
✗ Branch 53 not taken.
✗ Branch 54 not taken.
✓ Branch 55 taken 2 times.
✗ Branch 56 not taken.
✗ Branch 57 not taken.
✗ Branch 58 not taken.
✓ Branch 59 taken 2 times.
✗ Branch 60 not taken.
✗ Branch 61 not taken.
|
686 | PetscValidLogicalCollectiveInt(nep,maxits,3); |
| 276 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
|
686 | if (tol == (PetscReal)PETSC_DETERMINE) { |
| 277 | ✗ | nep->tol = PETSC_DETERMINE; | |
| 278 | ✗ | nep->state = NEP_STATE_INITIAL; | |
| 279 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
686 | } else if (tol != (PetscReal)PETSC_CURRENT) { |
| 280 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
676 | PetscCheck(tol>0.0,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of tol. Must be > 0"); |
| 281 | 676 | nep->tol = tol; | |
| 282 | } | ||
| 283 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
686 | if (maxits == PETSC_DETERMINE) { |
| 284 | 160 | nep->max_it = PETSC_DETERMINE; | |
| 285 | 160 | nep->state = NEP_STATE_INITIAL; | |
| 286 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
|
526 | } else if (maxits == PETSC_UNLIMITED) { |
| 287 | ✗ | nep->max_it = PETSC_INT_MAX; | |
| 288 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
526 | } else if (maxits != PETSC_CURRENT) { |
| 289 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
26 | PetscCheck(maxits>0,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of maxits. Must be > 0"); |
| 290 | 26 | nep->max_it = maxits; | |
| 291 | } | ||
| 292 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
132 | PetscFunctionReturn(PETSC_SUCCESS); |
| 293 | } | ||
| 294 | |||
| 295 | /*@ | ||
| 296 | NEPGetDimensions - Gets the number of eigenvalues to compute | ||
| 297 | and the dimension of the subspace. | ||
| 298 | |||
| 299 | Not Collective | ||
| 300 | |||
| 301 | Input Parameter: | ||
| 302 | . nep - the nonlinear eigensolver context | ||
| 303 | |||
| 304 | Output Parameters: | ||
| 305 | + nev - number of eigenvalues to compute | ||
| 306 | . ncv - the maximum dimension of the subspace to be used by the solver | ||
| 307 | - mpd - the maximum dimension allowed for the projected problem | ||
| 308 | |||
| 309 | Note: | ||
| 310 | The user can specify `NULL` for any parameter that is not needed. | ||
| 311 | |||
| 312 | Level: intermediate | ||
| 313 | |||
| 314 | .seealso: [](ch:nep), `NEPSetDimensions()` | ||
| 315 | @*/ | ||
| 316 | 539 | PetscErrorCode NEPGetDimensions(NEP nep,PetscInt *nev,PetscInt *ncv,PetscInt *mpd) | |
| 317 | { | ||
| 318 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
539 | PetscFunctionBegin; |
| 319 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
539 | PetscValidHeaderSpecific(nep,NEP_CLASSID,1); |
| 320 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
539 | if (nev) *nev = nep->nev; |
| 321 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
539 | if (ncv) *ncv = nep->ncv; |
| 322 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
539 | if (mpd) *mpd = nep->mpd; |
| 323 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
539 | PetscFunctionReturn(PETSC_SUCCESS); |
| 324 | } | ||
| 325 | |||
| 326 | /*@ | ||
| 327 | NEPSetDimensions - Sets the number of eigenvalues to compute | ||
| 328 | and the dimension of the subspace. | ||
| 329 | |||
| 330 | Logically Collective | ||
| 331 | |||
| 332 | Input Parameters: | ||
| 333 | + nep - the nonlinear eigensolver context | ||
| 334 | . nev - number of eigenvalues to compute | ||
| 335 | . ncv - the maximum dimension of the subspace to be used by the solver | ||
| 336 | - mpd - the maximum dimension allowed for the projected problem | ||
| 337 | |||
| 338 | Options Database Keys: | ||
| 339 | + -nep_nev \<nev\> - sets the number of eigenvalues | ||
| 340 | . -nep_ncv \<ncv\> - sets the dimension of the subspace | ||
| 341 | - -nep_mpd \<mpd\> - sets the maximum projected dimension | ||
| 342 | |||
| 343 | Notes: | ||
| 344 | Use `PETSC_DETERMINE` for `ncv` and `mpd` to assign a reasonably good value, which is | ||
| 345 | dependent on the solution method. For any of the arguments, use `PETSC_CURRENT` | ||
| 346 | to preserve the current value. | ||
| 347 | |||
| 348 | The parameters `ncv` and `mpd` are intimately related, so that the user is advised | ||
| 349 | to set one of them at most. Normal usage is\: | ||
| 350 | |||
| 351 | 1. in cases where `nev` is small, the user sets `ncv` (a reasonable default is `2*nev`). | ||
| 352 | 2. in cases where `nev` is large, the user sets `mpd`. | ||
| 353 | |||
| 354 | The value of `ncv` should always be between `nev` and `(nev+mpd)`, typically | ||
| 355 | `ncv=nev+mpd`. If `nev` is not too large, `mpd=nev` is a reasonable choice, otherwise | ||
| 356 | a smaller value should be used. | ||
| 357 | |||
| 358 | Level: intermediate | ||
| 359 | |||
| 360 | .seealso: [](ch:nep), `NEPGetDimensions()` | ||
| 361 | @*/ | ||
| 362 | 1034 | PetscErrorCode NEPSetDimensions(NEP nep,PetscInt nev,PetscInt ncv,PetscInt mpd) | |
| 363 | { | ||
| 364 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
1034 | PetscFunctionBegin; |
| 365 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
1034 | PetscValidHeaderSpecific(nep,NEP_CLASSID,1); |
| 366 |
27/62✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✓ Branch 13 taken 2 times.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
✓ Branch 16 taken 2 times.
✗ Branch 17 not taken.
✗ Branch 18 not taken.
✓ Branch 19 taken 2 times.
✓ Branch 20 taken 2 times.
✗ Branch 21 not taken.
✗ Branch 22 not taken.
✓ Branch 23 taken 2 times.
✓ Branch 24 taken 2 times.
✗ Branch 25 not taken.
✗ Branch 26 not taken.
✓ Branch 27 taken 2 times.
✓ Branch 28 taken 2 times.
✗ Branch 29 not taken.
✗ Branch 30 not taken.
✓ Branch 31 taken 2 times.
✓ Branch 32 taken 2 times.
✗ Branch 33 not taken.
✓ Branch 34 taken 2 times.
✗ Branch 35 not taken.
✗ Branch 36 not taken.
✓ Branch 37 taken 2 times.
✗ Branch 38 not taken.
✗ Branch 39 not taken.
✗ Branch 40 not taken.
✓ Branch 41 taken 2 times.
✗ Branch 42 not taken.
✗ Branch 43 not taken.
✗ Branch 44 not taken.
✓ Branch 45 taken 2 times.
✓ Branch 46 taken 2 times.
✗ Branch 47 not taken.
✗ Branch 48 not taken.
✓ Branch 49 taken 2 times.
✓ Branch 50 taken 2 times.
✗ Branch 51 not taken.
✓ Branch 52 taken 2 times.
✗ Branch 53 not taken.
✗ Branch 54 not taken.
✓ Branch 55 taken 2 times.
✗ Branch 56 not taken.
✗ Branch 57 not taken.
✗ Branch 58 not taken.
✓ Branch 59 taken 2 times.
✗ Branch 60 not taken.
✗ Branch 61 not taken.
|
1034 | PetscValidLogicalCollectiveInt(nep,nev,2); |
| 367 |
27/62✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✓ Branch 13 taken 2 times.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
✓ Branch 16 taken 2 times.
✗ Branch 17 not taken.
✗ Branch 18 not taken.
✓ Branch 19 taken 2 times.
✓ Branch 20 taken 2 times.
✗ Branch 21 not taken.
✗ Branch 22 not taken.
✓ Branch 23 taken 2 times.
✓ Branch 24 taken 2 times.
✗ Branch 25 not taken.
✗ Branch 26 not taken.
✓ Branch 27 taken 2 times.
✓ Branch 28 taken 2 times.
✗ Branch 29 not taken.
✗ Branch 30 not taken.
✓ Branch 31 taken 2 times.
✓ Branch 32 taken 2 times.
✗ Branch 33 not taken.
✓ Branch 34 taken 2 times.
✗ Branch 35 not taken.
✗ Branch 36 not taken.
✓ Branch 37 taken 2 times.
✗ Branch 38 not taken.
✗ Branch 39 not taken.
✗ Branch 40 not taken.
✓ Branch 41 taken 2 times.
✗ Branch 42 not taken.
✗ Branch 43 not taken.
✗ Branch 44 not taken.
✓ Branch 45 taken 2 times.
✓ Branch 46 taken 2 times.
✗ Branch 47 not taken.
✗ Branch 48 not taken.
✓ Branch 49 taken 2 times.
✓ Branch 50 taken 2 times.
✗ Branch 51 not taken.
✓ Branch 52 taken 2 times.
✗ Branch 53 not taken.
✗ Branch 54 not taken.
✓ Branch 55 taken 2 times.
✗ Branch 56 not taken.
✗ Branch 57 not taken.
✗ Branch 58 not taken.
✓ Branch 59 taken 2 times.
✗ Branch 60 not taken.
✗ Branch 61 not taken.
|
1034 | PetscValidLogicalCollectiveInt(nep,ncv,3); |
| 368 |
27/62✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✓ Branch 13 taken 2 times.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
✓ Branch 16 taken 2 times.
✗ Branch 17 not taken.
✗ Branch 18 not taken.
✓ Branch 19 taken 2 times.
✓ Branch 20 taken 2 times.
✗ Branch 21 not taken.
✗ Branch 22 not taken.
✓ Branch 23 taken 2 times.
✓ Branch 24 taken 2 times.
✗ Branch 25 not taken.
✗ Branch 26 not taken.
✓ Branch 27 taken 2 times.
✓ Branch 28 taken 2 times.
✗ Branch 29 not taken.
✗ Branch 30 not taken.
✓ Branch 31 taken 2 times.
✓ Branch 32 taken 2 times.
✗ Branch 33 not taken.
✓ Branch 34 taken 2 times.
✗ Branch 35 not taken.
✗ Branch 36 not taken.
✓ Branch 37 taken 2 times.
✗ Branch 38 not taken.
✗ Branch 39 not taken.
✗ Branch 40 not taken.
✓ Branch 41 taken 2 times.
✗ Branch 42 not taken.
✗ Branch 43 not taken.
✗ Branch 44 not taken.
✓ Branch 45 taken 2 times.
✓ Branch 46 taken 2 times.
✗ Branch 47 not taken.
✗ Branch 48 not taken.
✓ Branch 49 taken 2 times.
✓ Branch 50 taken 2 times.
✗ Branch 51 not taken.
✓ Branch 52 taken 2 times.
✗ Branch 53 not taken.
✗ Branch 54 not taken.
✓ Branch 55 taken 2 times.
✗ Branch 56 not taken.
✗ Branch 57 not taken.
✗ Branch 58 not taken.
✓ Branch 59 taken 2 times.
✗ Branch 60 not taken.
✗ Branch 61 not taken.
|
1034 | PetscValidLogicalCollectiveInt(nep,mpd,4); |
| 369 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
1034 | if (nev != PETSC_CURRENT) { |
| 370 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
1034 | PetscCheck(nev>0,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of nev. Must be > 0"); |
| 371 | 1034 | nep->nev = nev; | |
| 372 | } | ||
| 373 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
1034 | if (ncv == PETSC_DETERMINE) { |
| 374 | 878 | nep->ncv = PETSC_DETERMINE; | |
| 375 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 2 times.
|
156 | } else if (ncv != PETSC_CURRENT) { |
| 376 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
151 | PetscCheck(ncv>0,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of ncv. Must be > 0"); |
| 377 | 151 | nep->ncv = ncv; | |
| 378 | } | ||
| 379 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
1034 | if (mpd == PETSC_DETERMINE) { |
| 380 | 1009 | nep->mpd = PETSC_DETERMINE; | |
| 381 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
25 | } else if (mpd != PETSC_CURRENT) { |
| 382 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
10 | PetscCheck(mpd>0,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of mpd. Must be > 0"); |
| 383 | 10 | nep->mpd = mpd; | |
| 384 | } | ||
| 385 | 1034 | nep->state = NEP_STATE_INITIAL; | |
| 386 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
1034 | PetscFunctionReturn(PETSC_SUCCESS); |
| 387 | } | ||
| 388 | |||
| 389 | /*@ | ||
| 390 | NEPSetWhichEigenpairs - Specifies which portion of the spectrum is | ||
| 391 | to be sought. | ||
| 392 | |||
| 393 | Logically Collective | ||
| 394 | |||
| 395 | Input Parameters: | ||
| 396 | + nep - the nonlinear eigensolver context | ||
| 397 | - which - the portion of the spectrum to be sought, see `NEPWhich` for possible values | ||
| 398 | |||
| 399 | Options Database Keys: | ||
| 400 | + -nep_largest_magnitude - sets largest eigenvalues in magnitude | ||
| 401 | . -nep_smallest_magnitude - sets smallest eigenvalues in magnitude | ||
| 402 | . -nep_largest_real - sets largest real parts | ||
| 403 | . -nep_smallest_real - sets smallest real parts | ||
| 404 | . -nep_largest_imaginary - sets largest imaginary parts | ||
| 405 | . -nep_smallest_imaginary - sets smallest imaginary parts | ||
| 406 | . -nep_target_magnitude - sets eigenvalues closest to target | ||
| 407 | . -nep_target_real - sets real parts closest to target | ||
| 408 | . -nep_target_imaginary - sets imaginary parts closest to target | ||
| 409 | . -nep_all - sets all eigenvalues in a region | ||
| 410 | - -nep_which_user - select the user-defined selection criterion | ||
| 411 | |||
| 412 | Notes: | ||
| 413 | Not all eigensolvers implemented in `NEP` account for all the possible values | ||
| 414 | of `which`. Also, some values make sense only for certain types of | ||
| 415 | problems. If SLEPc is compiled for real numbers `NEP_LARGEST_IMAGINARY` | ||
| 416 | and `NEP_SMALLEST_IMAGINARY` use the absolute value of the imaginary part | ||
| 417 | for eigenvalue selection. | ||
| 418 | |||
| 419 | The target is a scalar value provided with `NEPSetTarget()`. | ||
| 420 | |||
| 421 | The criterion `NEP_TARGET_IMAGINARY` is available only in case PETSc and | ||
| 422 | SLEPc have been built with complex scalars. | ||
| 423 | |||
| 424 | `NEP_ALL` is intended for use in the context of the `NEPCISS` solver for | ||
| 425 | computing all eigenvalues in a region. | ||
| 426 | |||
| 427 | Level: intermediate | ||
| 428 | |||
| 429 | .seealso: [](ch:nep), `NEPGetWhichEigenpairs()`, `NEPSetTarget()`, `NEPSetDimensions()`, `NEPSetEigenvalueComparison()`, `NEPWhich` | ||
| 430 | @*/ | ||
| 431 | 791 | PetscErrorCode NEPSetWhichEigenpairs(NEP nep,NEPWhich which) | |
| 432 | { | ||
| 433 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
791 | PetscFunctionBegin; |
| 434 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
791 | PetscValidHeaderSpecific(nep,NEP_CLASSID,1); |
| 435 |
27/62✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✓ Branch 13 taken 2 times.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
✓ Branch 16 taken 2 times.
✗ Branch 17 not taken.
✗ Branch 18 not taken.
✓ Branch 19 taken 2 times.
✓ Branch 20 taken 2 times.
✗ Branch 21 not taken.
✗ Branch 22 not taken.
✓ Branch 23 taken 2 times.
✓ Branch 24 taken 2 times.
✗ Branch 25 not taken.
✗ Branch 26 not taken.
✓ Branch 27 taken 2 times.
✓ Branch 28 taken 2 times.
✗ Branch 29 not taken.
✗ Branch 30 not taken.
✓ Branch 31 taken 2 times.
✓ Branch 32 taken 2 times.
✗ Branch 33 not taken.
✓ Branch 34 taken 2 times.
✗ Branch 35 not taken.
✗ Branch 36 not taken.
✓ Branch 37 taken 2 times.
✗ Branch 38 not taken.
✗ Branch 39 not taken.
✗ Branch 40 not taken.
✓ Branch 41 taken 2 times.
✗ Branch 42 not taken.
✗ Branch 43 not taken.
✗ Branch 44 not taken.
✓ Branch 45 taken 2 times.
✓ Branch 46 taken 2 times.
✗ Branch 47 not taken.
✗ Branch 48 not taken.
✓ Branch 49 taken 2 times.
✓ Branch 50 taken 2 times.
✗ Branch 51 not taken.
✓ Branch 52 taken 2 times.
✗ Branch 53 not taken.
✗ Branch 54 not taken.
✓ Branch 55 taken 2 times.
✗ Branch 56 not taken.
✗ Branch 57 not taken.
✗ Branch 58 not taken.
✓ Branch 59 taken 2 times.
✗ Branch 60 not taken.
✗ Branch 61 not taken.
|
791 | PetscValidLogicalCollectiveEnum(nep,which,2); |
| 436 |
1/3✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
|
791 | switch (which) { |
| 437 | 791 | case NEP_LARGEST_MAGNITUDE: | |
| 438 | case NEP_SMALLEST_MAGNITUDE: | ||
| 439 | case NEP_LARGEST_REAL: | ||
| 440 | case NEP_SMALLEST_REAL: | ||
| 441 | case NEP_LARGEST_IMAGINARY: | ||
| 442 | case NEP_SMALLEST_IMAGINARY: | ||
| 443 | case NEP_TARGET_MAGNITUDE: | ||
| 444 | case NEP_TARGET_REAL: | ||
| 445 | #if defined(PETSC_USE_COMPLEX) | ||
| 446 | case NEP_TARGET_IMAGINARY: | ||
| 447 | #endif | ||
| 448 | case NEP_ALL: | ||
| 449 | case NEP_WHICH_USER: | ||
| 450 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
791 | if (nep->which != which) { |
| 451 | 771 | nep->state = NEP_STATE_INITIAL; | |
| 452 | 771 | nep->which = which; | |
| 453 | } | ||
| 454 | 152 | break; | |
| 455 | #if !defined(PETSC_USE_COMPLEX) | ||
| 456 | ✗ | case NEP_TARGET_IMAGINARY: | |
| 457 | ✗ | SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"NEP_TARGET_IMAGINARY can be used only with complex scalars"); | |
| 458 | #endif | ||
| 459 | ✗ | default: | |
| 460 | ✗ | SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'which' value"); | |
| 461 | } | ||
| 462 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
152 | PetscFunctionReturn(PETSC_SUCCESS); |
| 463 | } | ||
| 464 | |||
| 465 | /*@ | ||
| 466 | NEPGetWhichEigenpairs - Returns which portion of the spectrum is to be | ||
| 467 | sought. | ||
| 468 | |||
| 469 | Not Collective | ||
| 470 | |||
| 471 | Input Parameter: | ||
| 472 | . nep - the nonlinear eigensolver context | ||
| 473 | |||
| 474 | Output Parameter: | ||
| 475 | . which - the portion of the spectrum to be sought | ||
| 476 | |||
| 477 | Level: intermediate | ||
| 478 | |||
| 479 | .seealso: [](ch:nep), `NEPSetWhichEigenpairs()`, `NEPWhich` | ||
| 480 | @*/ | ||
| 481 | 16 | PetscErrorCode NEPGetWhichEigenpairs(NEP nep,NEPWhich *which) | |
| 482 | { | ||
| 483 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
16 | PetscFunctionBegin; |
| 484 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
16 | PetscValidHeaderSpecific(nep,NEP_CLASSID,1); |
| 485 |
2/8✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
16 | PetscAssertPointer(which,2); |
| 486 | 16 | *which = nep->which; | |
| 487 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
16 | PetscFunctionReturn(PETSC_SUCCESS); |
| 488 | } | ||
| 489 | |||
| 490 | /*@C | ||
| 491 | NEPSetEigenvalueComparison - Specifies the eigenvalue comparison function | ||
| 492 | when `NEPSetWhichEigenpairs()` is set to `NEP_WHICH_USER`. | ||
| 493 | |||
| 494 | Logically Collective | ||
| 495 | |||
| 496 | Input Parameters: | ||
| 497 | + nep - the nonlinear eigensolver context | ||
| 498 | . comp - a pointer to the comparison function, see `SlepcEigenvalueComparisonFn` for the calling sequence | ||
| 499 | - ctx - a context pointer (the last parameter to the comparison function) | ||
| 500 | |||
| 501 | Level: advanced | ||
| 502 | |||
| 503 | .seealso: [](ch:nep), `NEPSetWhichEigenpairs()`, `NEPWhich` | ||
| 504 | @*/ | ||
| 505 | 5 | PetscErrorCode NEPSetEigenvalueComparison(NEP nep,SlepcEigenvalueComparisonFn *comp,void *ctx) | |
| 506 | { | ||
| 507 |
1/2✓ Branch 0 taken 1 times.
✗ Branch 1 not taken.
|
5 | PetscFunctionBegin; |
| 508 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 1 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
5 | PetscValidHeaderSpecific(nep,NEP_CLASSID,1); |
| 509 | 5 | nep->sc->comparison = comp; | |
| 510 | 5 | nep->sc->comparisonctx = ctx; | |
| 511 | 5 | nep->which = NEP_WHICH_USER; | |
| 512 |
6/12✓ Branch 0 taken 1 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 1 times.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 1 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 1 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 1 times.
|
5 | PetscFunctionReturn(PETSC_SUCCESS); |
| 513 | } | ||
| 514 | |||
| 515 | /*@ | ||
| 516 | NEPSetProblemType - Specifies the type of the nonlinear eigenvalue problem. | ||
| 517 | |||
| 518 | Logically Collective | ||
| 519 | |||
| 520 | Input Parameters: | ||
| 521 | + nep - the nonlinear eigensolver context | ||
| 522 | - type - a known type of nonlinear eigenvalue problem | ||
| 523 | |||
| 524 | Options Database Keys: | ||
| 525 | + -nep_general - general problem with no particular structure | ||
| 526 | - -nep_rational - a rational eigenvalue problem defined in split form with all $f_i$ rational | ||
| 527 | |||
| 528 | Notes: | ||
| 529 | See `NEPProblemType` for possible problem types. | ||
| 530 | |||
| 531 | This function is used to provide a hint to the `NEP` solver to exploit certain | ||
| 532 | properties of the nonlinear eigenproblem. This hint may be used or not, | ||
| 533 | depending on the solver. By default, no particular structure is assumed. | ||
| 534 | |||
| 535 | Level: intermediate | ||
| 536 | |||
| 537 | .seealso: [](ch:nep), `NEPSetType()`, `NEPGetProblemType()`, `NEPProblemType` | ||
| 538 | @*/ | ||
| 539 | 1305 | PetscErrorCode NEPSetProblemType(NEP nep,NEPProblemType type) | |
| 540 | { | ||
| 541 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
1305 | PetscFunctionBegin; |
| 542 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
1305 | PetscValidHeaderSpecific(nep,NEP_CLASSID,1); |
| 543 |
27/62✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✓ Branch 13 taken 2 times.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
✓ Branch 16 taken 2 times.
✗ Branch 17 not taken.
✗ Branch 18 not taken.
✓ Branch 19 taken 2 times.
✓ Branch 20 taken 2 times.
✗ Branch 21 not taken.
✗ Branch 22 not taken.
✓ Branch 23 taken 2 times.
✓ Branch 24 taken 2 times.
✗ Branch 25 not taken.
✗ Branch 26 not taken.
✓ Branch 27 taken 2 times.
✓ Branch 28 taken 2 times.
✗ Branch 29 not taken.
✗ Branch 30 not taken.
✓ Branch 31 taken 2 times.
✓ Branch 32 taken 2 times.
✗ Branch 33 not taken.
✓ Branch 34 taken 2 times.
✗ Branch 35 not taken.
✗ Branch 36 not taken.
✓ Branch 37 taken 2 times.
✗ Branch 38 not taken.
✗ Branch 39 not taken.
✗ Branch 40 not taken.
✓ Branch 41 taken 2 times.
✗ Branch 42 not taken.
✗ Branch 43 not taken.
✗ Branch 44 not taken.
✓ Branch 45 taken 2 times.
✓ Branch 46 taken 2 times.
✗ Branch 47 not taken.
✗ Branch 48 not taken.
✓ Branch 49 taken 2 times.
✓ Branch 50 taken 2 times.
✗ Branch 51 not taken.
✓ Branch 52 taken 2 times.
✗ Branch 53 not taken.
✗ Branch 54 not taken.
✓ Branch 55 taken 2 times.
✗ Branch 56 not taken.
✗ Branch 57 not taken.
✗ Branch 58 not taken.
✓ Branch 59 taken 2 times.
✗ Branch 60 not taken.
✗ Branch 61 not taken.
|
1305 | PetscValidLogicalCollectiveEnum(nep,type,2); |
| 544 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
1305 | PetscCheck(type==NEP_GENERAL || type==NEP_RATIONAL,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_WRONG,"Unknown eigenvalue problem type"); |
| 545 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
1305 | if (type != nep->problem_type) { |
| 546 | 1305 | nep->problem_type = type; | |
| 547 | 1305 | nep->state = NEP_STATE_INITIAL; | |
| 548 | } | ||
| 549 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
247 | PetscFunctionReturn(PETSC_SUCCESS); |
| 550 | } | ||
| 551 | |||
| 552 | /*@ | ||
| 553 | NEPGetProblemType - Gets the problem type from the `NEP` object. | ||
| 554 | |||
| 555 | Not Collective | ||
| 556 | |||
| 557 | Input Parameter: | ||
| 558 | . nep - the nonlinear eigensolver context | ||
| 559 | |||
| 560 | Output Parameter: | ||
| 561 | . type - the problem type | ||
| 562 | |||
| 563 | Level: intermediate | ||
| 564 | |||
| 565 | .seealso: [](ch:nep), `NEPSetProblemType()`, `NEPProblemType` | ||
| 566 | @*/ | ||
| 567 | 32 | PetscErrorCode NEPGetProblemType(NEP nep,NEPProblemType *type) | |
| 568 | { | ||
| 569 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
32 | PetscFunctionBegin; |
| 570 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
32 | PetscValidHeaderSpecific(nep,NEP_CLASSID,1); |
| 571 |
2/8✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
32 | PetscAssertPointer(type,2); |
| 572 | 32 | *type = nep->problem_type; | |
| 573 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
32 | PetscFunctionReturn(PETSC_SUCCESS); |
| 574 | } | ||
| 575 | |||
| 576 | /*@ | ||
| 577 | NEPSetTwoSided - Sets the solver to use a two-sided variant so that left | ||
| 578 | eigenvectors are also computed. | ||
| 579 | |||
| 580 | Logically Collective | ||
| 581 | |||
| 582 | Input Parameters: | ||
| 583 | + nep - the nonlinear eigensolver context | ||
| 584 | - twosided - whether the two-sided variant is to be used or not | ||
| 585 | |||
| 586 | Options Database Key: | ||
| 587 | . -nep_two_sided - toggles the twosided flag | ||
| 588 | |||
| 589 | Notes: | ||
| 590 | If the user sets `twosided`=`PETSC_TRUE` then the solver uses a variant of | ||
| 591 | the algorithm that computes both right and left eigenvectors. This is | ||
| 592 | usually much more costly. This option is not available in all solvers, | ||
| 593 | see table [](#tab:solversn). | ||
| 594 | |||
| 595 | When using two-sided solvers, the problem matrices must have both the | ||
| 596 | `MATOP_MULT` and `MATOP_MULT_TRANSPOSE` operations defined. | ||
| 597 | |||
| 598 | Level: advanced | ||
| 599 | |||
| 600 | .seealso: [](ch:nep), `NEPGetTwoSided()`, `NEPGetLeftEigenvector()` | ||
| 601 | @*/ | ||
| 602 | 100 | PetscErrorCode NEPSetTwoSided(NEP nep,PetscBool twosided) | |
| 603 | { | ||
| 604 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
100 | PetscFunctionBegin; |
| 605 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
100 | PetscValidHeaderSpecific(nep,NEP_CLASSID,1); |
| 606 |
27/62✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✓ Branch 13 taken 2 times.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
✓ Branch 16 taken 2 times.
✗ Branch 17 not taken.
✗ Branch 18 not taken.
✓ Branch 19 taken 2 times.
✓ Branch 20 taken 2 times.
✗ Branch 21 not taken.
✗ Branch 22 not taken.
✓ Branch 23 taken 2 times.
✓ Branch 24 taken 2 times.
✗ Branch 25 not taken.
✗ Branch 26 not taken.
✓ Branch 27 taken 2 times.
✓ Branch 28 taken 2 times.
✗ Branch 29 not taken.
✗ Branch 30 not taken.
✓ Branch 31 taken 2 times.
✓ Branch 32 taken 2 times.
✗ Branch 33 not taken.
✓ Branch 34 taken 2 times.
✗ Branch 35 not taken.
✗ Branch 36 not taken.
✓ Branch 37 taken 2 times.
✗ Branch 38 not taken.
✗ Branch 39 not taken.
✗ Branch 40 not taken.
✓ Branch 41 taken 2 times.
✗ Branch 42 not taken.
✗ Branch 43 not taken.
✗ Branch 44 not taken.
✓ Branch 45 taken 2 times.
✓ Branch 46 taken 2 times.
✗ Branch 47 not taken.
✗ Branch 48 not taken.
✓ Branch 49 taken 2 times.
✓ Branch 50 taken 2 times.
✗ Branch 51 not taken.
✓ Branch 52 taken 2 times.
✗ Branch 53 not taken.
✗ Branch 54 not taken.
✓ Branch 55 taken 2 times.
✗ Branch 56 not taken.
✗ Branch 57 not taken.
✗ Branch 58 not taken.
✓ Branch 59 taken 2 times.
✗ Branch 60 not taken.
✗ Branch 61 not taken.
|
100 | PetscValidLogicalCollectiveBool(nep,twosided,2); |
| 607 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
100 | if (twosided!=nep->twosided) { |
| 608 | 70 | nep->twosided = twosided; | |
| 609 | 70 | nep->state = NEP_STATE_INITIAL; | |
| 610 | } | ||
| 611 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
100 | PetscFunctionReturn(PETSC_SUCCESS); |
| 612 | } | ||
| 613 | |||
| 614 | /*@ | ||
| 615 | NEPGetTwoSided - Returns the flag indicating whether a two-sided variant | ||
| 616 | of the algorithm is being used or not. | ||
| 617 | |||
| 618 | Not Collective | ||
| 619 | |||
| 620 | Input Parameter: | ||
| 621 | . nep - the nonlinear eigensolver context | ||
| 622 | |||
| 623 | Output Parameter: | ||
| 624 | . twosided - the returned flag | ||
| 625 | |||
| 626 | Level: advanced | ||
| 627 | |||
| 628 | .seealso: [](ch:nep), `NEPSetTwoSided()` | ||
| 629 | @*/ | ||
| 630 | 30 | PetscErrorCode NEPGetTwoSided(NEP nep,PetscBool *twosided) | |
| 631 | { | ||
| 632 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
30 | PetscFunctionBegin; |
| 633 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
30 | PetscValidHeaderSpecific(nep,NEP_CLASSID,1); |
| 634 |
2/8✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
30 | PetscAssertPointer(twosided,2); |
| 635 | 30 | *twosided = nep->twosided; | |
| 636 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
30 | PetscFunctionReturn(PETSC_SUCCESS); |
| 637 | } | ||
| 638 | |||
| 639 | /*@C | ||
| 640 | NEPSetConvergenceTestFunction - Sets a function to compute the error estimate | ||
| 641 | used in the convergence test. | ||
| 642 | |||
| 643 | Logically Collective | ||
| 644 | |||
| 645 | Input Parameters: | ||
| 646 | + nep - the nonlinear eigensolver context | ||
| 647 | . conv - convergence test function, see `NEPConvergenceTestFn` for the calling sequence | ||
| 648 | . ctx - context for private data for the convergence routine (may be `NULL`) | ||
| 649 | - destroy - a routine for destroying the context (may be `NULL`), see `PetscCtxDestroyFn` | ||
| 650 | for the calling sequence | ||
| 651 | |||
| 652 | Notes: | ||
| 653 | When this is called with a user-defined function, then the convergence | ||
| 654 | criterion is set to `NEP_CONV_USER`, see `NEPSetConvergenceTest()`. | ||
| 655 | |||
| 656 | If the error estimate returned by the convergence test function is less than | ||
| 657 | the tolerance, then the eigenvalue is accepted as converged. | ||
| 658 | |||
| 659 | Level: advanced | ||
| 660 | |||
| 661 | .seealso: [](ch:nep), `NEPSetConvergenceTest()`, `NEPSetTolerances()` | ||
| 662 | @*/ | ||
| 663 | 10 | PetscErrorCode NEPSetConvergenceTestFunction(NEP nep,NEPConvergenceTestFn *conv,void *ctx,PetscCtxDestroyFn *destroy) | |
| 664 | { | ||
| 665 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
10 | PetscFunctionBegin; |
| 666 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
10 | PetscValidHeaderSpecific(nep,NEP_CLASSID,1); |
| 667 |
1/8✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
10 | if (nep->convergeddestroy) PetscCall((*nep->convergeddestroy)(&nep->convergedctx)); |
| 668 | 10 | nep->convergeduser = conv; | |
| 669 | 10 | nep->convergeddestroy = destroy; | |
| 670 | 10 | nep->convergedctx = ctx; | |
| 671 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
|
10 | if (conv == NEPConvergedRelative) nep->conv = NEP_CONV_REL; |
| 672 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
|
10 | else if (conv == NEPConvergedNorm) nep->conv = NEP_CONV_NORM; |
| 673 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
|
10 | else if (conv == NEPConvergedAbsolute) nep->conv = NEP_CONV_ABS; |
| 674 | else { | ||
| 675 | 10 | nep->conv = NEP_CONV_USER; | |
| 676 | 10 | nep->converged = nep->convergeduser; | |
| 677 | } | ||
| 678 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
2 | PetscFunctionReturn(PETSC_SUCCESS); |
| 679 | } | ||
| 680 | |||
| 681 | /*@ | ||
| 682 | NEPSetConvergenceTest - Specifies how to compute the error estimate | ||
| 683 | used in the convergence test. | ||
| 684 | |||
| 685 | Logically Collective | ||
| 686 | |||
| 687 | Input Parameters: | ||
| 688 | + nep - the nonlinear eigensolver context | ||
| 689 | - conv - the type of convergence test, see `NEPConv` for possible values | ||
| 690 | |||
| 691 | Options Database Keys: | ||
| 692 | + -nep_conv_abs - sets the absolute convergence test | ||
| 693 | . -nep_conv_rel - sets the convergence test relative to the eigenvalue | ||
| 694 | . -nep_conv_norm - sets the convergence test relative to the matrix norms | ||
| 695 | - -nep_conv_user - selects the user-defined convergence test | ||
| 696 | |||
| 697 | Level: intermediate | ||
| 698 | |||
| 699 | .seealso: [](ch:nep), `NEPGetConvergenceTest()`, `NEPSetConvergenceTestFunction()`, `NEPSetStoppingTest()`, `NEPConv` | ||
| 700 | @*/ | ||
| 701 | 26 | PetscErrorCode NEPSetConvergenceTest(NEP nep,NEPConv conv) | |
| 702 | { | ||
| 703 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
26 | PetscFunctionBegin; |
| 704 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
26 | PetscValidHeaderSpecific(nep,NEP_CLASSID,1); |
| 705 |
27/62✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✓ Branch 13 taken 2 times.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
✓ Branch 16 taken 2 times.
✗ Branch 17 not taken.
✗ Branch 18 not taken.
✓ Branch 19 taken 2 times.
✓ Branch 20 taken 2 times.
✗ Branch 21 not taken.
✗ Branch 22 not taken.
✓ Branch 23 taken 2 times.
✓ Branch 24 taken 2 times.
✗ Branch 25 not taken.
✗ Branch 26 not taken.
✓ Branch 27 taken 2 times.
✓ Branch 28 taken 2 times.
✗ Branch 29 not taken.
✗ Branch 30 not taken.
✓ Branch 31 taken 2 times.
✓ Branch 32 taken 2 times.
✗ Branch 33 not taken.
✓ Branch 34 taken 2 times.
✗ Branch 35 not taken.
✗ Branch 36 not taken.
✓ Branch 37 taken 2 times.
✗ Branch 38 not taken.
✗ Branch 39 not taken.
✗ Branch 40 not taken.
✓ Branch 41 taken 2 times.
✗ Branch 42 not taken.
✗ Branch 43 not taken.
✗ Branch 44 not taken.
✓ Branch 45 taken 2 times.
✓ Branch 46 taken 2 times.
✗ Branch 47 not taken.
✗ Branch 48 not taken.
✓ Branch 49 taken 2 times.
✓ Branch 50 taken 2 times.
✗ Branch 51 not taken.
✓ Branch 52 taken 2 times.
✗ Branch 53 not taken.
✗ Branch 54 not taken.
✓ Branch 55 taken 2 times.
✗ Branch 56 not taken.
✗ Branch 57 not taken.
✗ Branch 58 not taken.
✓ Branch 59 taken 2 times.
✗ Branch 60 not taken.
✗ Branch 61 not taken.
|
26 | PetscValidLogicalCollectiveEnum(nep,conv,2); |
| 706 |
2/5✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 10 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
|
26 | switch (conv) { |
| 707 | 16 | case NEP_CONV_ABS: nep->converged = NEPConvergedAbsolute; break; | |
| 708 | ✗ | case NEP_CONV_REL: nep->converged = NEPConvergedRelative; break; | |
| 709 | 10 | case NEP_CONV_NORM: nep->converged = NEPConvergedNorm; break; | |
| 710 | ✗ | case NEP_CONV_USER: | |
| 711 | ✗ | PetscCheck(nep->convergeduser,PetscObjectComm((PetscObject)nep),PETSC_ERR_ORDER,"Must call NEPSetConvergenceTestFunction() first"); | |
| 712 | ✗ | nep->converged = nep->convergeduser; | |
| 713 | ✗ | break; | |
| 714 | ✗ | default: | |
| 715 | ✗ | SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'conv' value"); | |
| 716 | } | ||
| 717 | 26 | nep->conv = conv; | |
| 718 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
26 | PetscFunctionReturn(PETSC_SUCCESS); |
| 719 | } | ||
| 720 | |||
| 721 | /*@ | ||
| 722 | NEPGetConvergenceTest - Gets the method used to compute the error estimate | ||
| 723 | used in the convergence test. | ||
| 724 | |||
| 725 | Not Collective | ||
| 726 | |||
| 727 | Input Parameter: | ||
| 728 | . nep - the nonlinear eigensolver context | ||
| 729 | |||
| 730 | Output Parameter: | ||
| 731 | . conv - the type of convergence test | ||
| 732 | |||
| 733 | Level: intermediate | ||
| 734 | |||
| 735 | .seealso: [](ch:nep), `NEPSetConvergenceTest()`, `NEPConv` | ||
| 736 | @*/ | ||
| 737 | 16 | PetscErrorCode NEPGetConvergenceTest(NEP nep,NEPConv *conv) | |
| 738 | { | ||
| 739 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
16 | PetscFunctionBegin; |
| 740 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
16 | PetscValidHeaderSpecific(nep,NEP_CLASSID,1); |
| 741 |
2/8✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
16 | PetscAssertPointer(conv,2); |
| 742 | 16 | *conv = nep->conv; | |
| 743 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
16 | PetscFunctionReturn(PETSC_SUCCESS); |
| 744 | } | ||
| 745 | |||
| 746 | /*@C | ||
| 747 | NEPSetStoppingTestFunction - Sets a function to decide when to stop the outer | ||
| 748 | iteration of the eigensolver. | ||
| 749 | |||
| 750 | Logically Collective | ||
| 751 | |||
| 752 | Input Parameters: | ||
| 753 | + nep - the nonlinear eigensolver context | ||
| 754 | . stop - the stopping test function, see `NEPStoppingTestFn` for the calling sequence | ||
| 755 | . ctx - context for private data for the stopping routine (may be `NULL`) | ||
| 756 | - destroy - a routine for destroying the context (may be `NULL`), see `PetscCtxDestroyFn` | ||
| 757 | for the calling sequence | ||
| 758 | |||
| 759 | Note: | ||
| 760 | When implementing a function for this, normal usage is to first call the | ||
| 761 | default routine `NEPStoppingBasic()` and then set `reason` to `NEP_CONVERGED_USER` | ||
| 762 | if some user-defined conditions have been met. To let the eigensolver continue | ||
| 763 | iterating, the result must be left as `NEP_CONVERGED_ITERATING`. | ||
| 764 | |||
| 765 | Level: advanced | ||
| 766 | |||
| 767 | .seealso: [](ch:nep), `NEPSetStoppingTest()`, `NEPStoppingBasic()` | ||
| 768 | @*/ | ||
| 769 | 10 | PetscErrorCode NEPSetStoppingTestFunction(NEP nep,NEPStoppingTestFn *stop,void *ctx,PetscCtxDestroyFn *destroy) | |
| 770 | { | ||
| 771 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
10 | PetscFunctionBegin; |
| 772 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
10 | PetscValidHeaderSpecific(nep,NEP_CLASSID,1); |
| 773 |
1/8✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
10 | if (nep->stoppingdestroy) PetscCall((*nep->stoppingdestroy)(&nep->stoppingctx)); |
| 774 | 10 | nep->stoppinguser = stop; | |
| 775 | 10 | nep->stoppingdestroy = destroy; | |
| 776 | 10 | nep->stoppingctx = ctx; | |
| 777 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
|
10 | if (stop == NEPStoppingBasic) nep->stop = NEP_STOP_BASIC; |
| 778 | else { | ||
| 779 | 10 | nep->stop = NEP_STOP_USER; | |
| 780 | 10 | nep->stopping = nep->stoppinguser; | |
| 781 | } | ||
| 782 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
2 | PetscFunctionReturn(PETSC_SUCCESS); |
| 783 | } | ||
| 784 | |||
| 785 | /*@ | ||
| 786 | NEPSetStoppingTest - Specifies how to decide the termination of the outer | ||
| 787 | loop of the eigensolver. | ||
| 788 | |||
| 789 | Logically Collective | ||
| 790 | |||
| 791 | Input Parameters: | ||
| 792 | + nep - the nonlinear eigensolver context | ||
| 793 | - stop - the type of stopping test, see `NEPStop` | ||
| 794 | |||
| 795 | Options Database Keys: | ||
| 796 | + -nep_stop_basic - sets the default stopping test | ||
| 797 | - -nep_stop_user - selects the user-defined stopping test | ||
| 798 | |||
| 799 | Level: advanced | ||
| 800 | |||
| 801 | .seealso: [](ch:nep), `NEPGetStoppingTest()`, `NEPSetStoppingTestFunction()`, `NEPSetConvergenceTest()`, `NEPStop` | ||
| 802 | @*/ | ||
| 803 | 16 | PetscErrorCode NEPSetStoppingTest(NEP nep,NEPStop stop) | |
| 804 | { | ||
| 805 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
16 | PetscFunctionBegin; |
| 806 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
16 | PetscValidHeaderSpecific(nep,NEP_CLASSID,1); |
| 807 |
27/62✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✓ Branch 13 taken 2 times.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
✓ Branch 16 taken 2 times.
✗ Branch 17 not taken.
✗ Branch 18 not taken.
✓ Branch 19 taken 2 times.
✓ Branch 20 taken 2 times.
✗ Branch 21 not taken.
✗ Branch 22 not taken.
✓ Branch 23 taken 2 times.
✓ Branch 24 taken 2 times.
✗ Branch 25 not taken.
✗ Branch 26 not taken.
✓ Branch 27 taken 2 times.
✓ Branch 28 taken 2 times.
✗ Branch 29 not taken.
✗ Branch 30 not taken.
✓ Branch 31 taken 2 times.
✓ Branch 32 taken 2 times.
✗ Branch 33 not taken.
✓ Branch 34 taken 2 times.
✗ Branch 35 not taken.
✗ Branch 36 not taken.
✓ Branch 37 taken 2 times.
✗ Branch 38 not taken.
✗ Branch 39 not taken.
✗ Branch 40 not taken.
✓ Branch 41 taken 2 times.
✗ Branch 42 not taken.
✗ Branch 43 not taken.
✗ Branch 44 not taken.
✓ Branch 45 taken 2 times.
✓ Branch 46 taken 2 times.
✗ Branch 47 not taken.
✗ Branch 48 not taken.
✓ Branch 49 taken 2 times.
✓ Branch 50 taken 2 times.
✗ Branch 51 not taken.
✓ Branch 52 taken 2 times.
✗ Branch 53 not taken.
✗ Branch 54 not taken.
✓ Branch 55 taken 2 times.
✗ Branch 56 not taken.
✗ Branch 57 not taken.
✗ Branch 58 not taken.
✓ Branch 59 taken 2 times.
✗ Branch 60 not taken.
✗ Branch 61 not taken.
|
16 | PetscValidLogicalCollectiveEnum(nep,stop,2); |
| 808 |
1/3✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
|
16 | switch (stop) { |
| 809 | 16 | case NEP_STOP_BASIC: nep->stopping = NEPStoppingBasic; break; | |
| 810 | ✗ | case NEP_STOP_USER: | |
| 811 | ✗ | PetscCheck(nep->stoppinguser,PetscObjectComm((PetscObject)nep),PETSC_ERR_ORDER,"Must call NEPSetStoppingTestFunction() first"); | |
| 812 | ✗ | nep->stopping = nep->stoppinguser; | |
| 813 | ✗ | break; | |
| 814 | ✗ | default: | |
| 815 | ✗ | SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'stop' value"); | |
| 816 | } | ||
| 817 | 16 | nep->stop = stop; | |
| 818 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
16 | PetscFunctionReturn(PETSC_SUCCESS); |
| 819 | } | ||
| 820 | |||
| 821 | /*@ | ||
| 822 | NEPGetStoppingTest - Gets the method used to decide the termination of the outer | ||
| 823 | loop of the eigensolver. | ||
| 824 | |||
| 825 | Not Collective | ||
| 826 | |||
| 827 | Input Parameter: | ||
| 828 | . nep - the nonlinear eigensolver context | ||
| 829 | |||
| 830 | Output Parameter: | ||
| 831 | . stop - the type of stopping test | ||
| 832 | |||
| 833 | Level: advanced | ||
| 834 | |||
| 835 | .seealso: [](ch:nep), `NEPSetStoppingTest()`, `NEPStop` | ||
| 836 | @*/ | ||
| 837 | 16 | PetscErrorCode NEPGetStoppingTest(NEP nep,NEPStop *stop) | |
| 838 | { | ||
| 839 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
16 | PetscFunctionBegin; |
| 840 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
16 | PetscValidHeaderSpecific(nep,NEP_CLASSID,1); |
| 841 |
2/8✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
16 | PetscAssertPointer(stop,2); |
| 842 | 16 | *stop = nep->stop; | |
| 843 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
16 | PetscFunctionReturn(PETSC_SUCCESS); |
| 844 | } | ||
| 845 | |||
| 846 | /*@ | ||
| 847 | NEPSetTrackAll - Specifies if the solver must compute the residual of all | ||
| 848 | approximate eigenpairs or not. | ||
| 849 | |||
| 850 | Logically Collective | ||
| 851 | |||
| 852 | Input Parameters: | ||
| 853 | + nep - the nonlinear eigensolver context | ||
| 854 | - trackall - whether compute all residuals or not | ||
| 855 | |||
| 856 | Notes: | ||
| 857 | If the user sets `trackall`=`PETSC_TRUE` then the solver explicitly computes | ||
| 858 | the residual for each eigenpair approximation. Computing the residual is | ||
| 859 | usually an expensive operation and solvers commonly compute the associated | ||
| 860 | residual to the first unconverged eigenpair. | ||
| 861 | |||
| 862 | The option `-nep_monitor_all` automatically activates this option. | ||
| 863 | |||
| 864 | Level: developer | ||
| 865 | |||
| 866 | .seealso: [](ch:nep), `NEPGetTrackAll()` | ||
| 867 | @*/ | ||
| 868 | 18 | PetscErrorCode NEPSetTrackAll(NEP nep,PetscBool trackall) | |
| 869 | { | ||
| 870 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
18 | PetscFunctionBegin; |
| 871 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
18 | PetscValidHeaderSpecific(nep,NEP_CLASSID,1); |
| 872 |
27/62✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✓ Branch 13 taken 2 times.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
✓ Branch 16 taken 2 times.
✗ Branch 17 not taken.
✗ Branch 18 not taken.
✓ Branch 19 taken 2 times.
✓ Branch 20 taken 2 times.
✗ Branch 21 not taken.
✗ Branch 22 not taken.
✓ Branch 23 taken 2 times.
✓ Branch 24 taken 2 times.
✗ Branch 25 not taken.
✗ Branch 26 not taken.
✓ Branch 27 taken 2 times.
✓ Branch 28 taken 2 times.
✗ Branch 29 not taken.
✗ Branch 30 not taken.
✓ Branch 31 taken 2 times.
✓ Branch 32 taken 2 times.
✗ Branch 33 not taken.
✓ Branch 34 taken 2 times.
✗ Branch 35 not taken.
✗ Branch 36 not taken.
✓ Branch 37 taken 2 times.
✗ Branch 38 not taken.
✗ Branch 39 not taken.
✗ Branch 40 not taken.
✓ Branch 41 taken 2 times.
✗ Branch 42 not taken.
✗ Branch 43 not taken.
✗ Branch 44 not taken.
✓ Branch 45 taken 2 times.
✓ Branch 46 taken 2 times.
✗ Branch 47 not taken.
✗ Branch 48 not taken.
✓ Branch 49 taken 2 times.
✓ Branch 50 taken 2 times.
✗ Branch 51 not taken.
✓ Branch 52 taken 2 times.
✗ Branch 53 not taken.
✗ Branch 54 not taken.
✓ Branch 55 taken 2 times.
✗ Branch 56 not taken.
✗ Branch 57 not taken.
✗ Branch 58 not taken.
✓ Branch 59 taken 2 times.
✗ Branch 60 not taken.
✗ Branch 61 not taken.
|
18 | PetscValidLogicalCollectiveBool(nep,trackall,2); |
| 873 | 18 | nep->trackall = trackall; | |
| 874 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
18 | PetscFunctionReturn(PETSC_SUCCESS); |
| 875 | } | ||
| 876 | |||
| 877 | /*@ | ||
| 878 | NEPGetTrackAll - Returns the flag indicating whether all residual norms must | ||
| 879 | be computed or not. | ||
| 880 | |||
| 881 | Not Collective | ||
| 882 | |||
| 883 | Input Parameter: | ||
| 884 | . nep - the nonlinear eigensolver context | ||
| 885 | |||
| 886 | Output Parameter: | ||
| 887 | . trackall - the returned flag | ||
| 888 | |||
| 889 | Level: developer | ||
| 890 | |||
| 891 | .seealso: [](ch:nep), `NEPSetTrackAll()` | ||
| 892 | @*/ | ||
| 893 | 280 | PetscErrorCode NEPGetTrackAll(NEP nep,PetscBool *trackall) | |
| 894 | { | ||
| 895 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
280 | PetscFunctionBegin; |
| 896 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
280 | PetscValidHeaderSpecific(nep,NEP_CLASSID,1); |
| 897 |
2/8✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
280 | PetscAssertPointer(trackall,2); |
| 898 | 280 | *trackall = nep->trackall; | |
| 899 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
280 | PetscFunctionReturn(PETSC_SUCCESS); |
| 900 | } | ||
| 901 | |||
| 902 | /*@ | ||
| 903 | NEPSetRefine - Specifies the refinement type (and options) to be used | ||
| 904 | after the solve. | ||
| 905 | |||
| 906 | Logically Collective | ||
| 907 | |||
| 908 | Input Parameters: | ||
| 909 | + nep - the nonlinear eigensolver context | ||
| 910 | . refine - refinement type, see `NEPRefine` for possible values | ||
| 911 | . npart - number of partitions of the communicator | ||
| 912 | . tol - the convergence tolerance | ||
| 913 | . its - maximum number of refinement iterations | ||
| 914 | - scheme - which scheme to be used for solving the involved linear systems, see `NEPRefineScheme` | ||
| 915 | for possible values | ||
| 916 | |||
| 917 | Options Database Keys: | ||
| 918 | + -nep_refine \<refine\> - set the refinement type, one of `none`,`simple`,`multiple` | ||
| 919 | . -nep_refine_partitions \<npart\> - set the number of partitions | ||
| 920 | . -nep_refine_tol \<tol\> - set the tolerance | ||
| 921 | . -nep_refine_its \<its\> - set the number of iterations | ||
| 922 | - -nep_refine_scheme \<scheme\> - set the scheme for the linear solves, `schur`,`mbe`, or `explicit` | ||
| 923 | |||
| 924 | Notes: | ||
| 925 | This function configures the parameters of Newton iterative refinement, | ||
| 926 | see section [](#sec:refine) for a discussion of the different strategies | ||
| 927 | in the context of polynomial eigenproblems. | ||
| 928 | |||
| 929 | By default, iterative refinement is disabled, since it may be very | ||
| 930 | costly. There are two possible refinement strategies, simple and multiple. | ||
| 931 | The simple approach performs iterative refinement on each of the | ||
| 932 | converged eigenpairs individually, whereas the multiple strategy works | ||
| 933 | with the invariant pair as a whole, refining all eigenpairs simultaneously. | ||
| 934 | The latter may be required for the case of multiple eigenvalues. | ||
| 935 | |||
| 936 | In some cases, especially when using direct solvers within the | ||
| 937 | iterative refinement method, it may be helpful for improved scalability | ||
| 938 | to split the communicator in several partitions. The `npart` parameter | ||
| 939 | indicates how many partitions to use (defaults to 1). | ||
| 940 | |||
| 941 | The `tol` and `its` parameters specify the stopping criterion. In the simple | ||
| 942 | method, refinement continues until the residual of each eigenpair is | ||
| 943 | below the tolerance (`tol` defaults to the `NEP` tolerance, but may be set to a | ||
| 944 | different value). In contrast, the multiple method simply performs its | ||
| 945 | refinement iterations (just one by default). | ||
| 946 | |||
| 947 | The `scheme` argument is used to change the way in which linear systems are | ||
| 948 | solved. Possible choices are explicit, mixed block elimination (MBE), | ||
| 949 | and Schur complement. | ||
| 950 | |||
| 951 | Use `PETSC_CURRENT` to retain the current value of `npart`, `tol` or `its`. Use | ||
| 952 | `PETSC_DETERMINE` to assign a default value. | ||
| 953 | |||
| 954 | Level: intermediate | ||
| 955 | |||
| 956 | .seealso: [](ch:nep), [](#sec:refine), `NEPGetRefine()` | ||
| 957 | @*/ | ||
| 958 | 136 | PetscErrorCode NEPSetRefine(NEP nep,NEPRefine refine,PetscInt npart,PetscReal tol,PetscInt its,NEPRefineScheme scheme) | |
| 959 | { | ||
| 960 | 136 | PetscMPIInt size; | |
| 961 | |||
| 962 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
136 | PetscFunctionBegin; |
| 963 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
136 | PetscValidHeaderSpecific(nep,NEP_CLASSID,1); |
| 964 |
27/62✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✓ Branch 13 taken 2 times.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
✓ Branch 16 taken 2 times.
✗ Branch 17 not taken.
✗ Branch 18 not taken.
✓ Branch 19 taken 2 times.
✓ Branch 20 taken 2 times.
✗ Branch 21 not taken.
✗ Branch 22 not taken.
✓ Branch 23 taken 2 times.
✓ Branch 24 taken 2 times.
✗ Branch 25 not taken.
✗ Branch 26 not taken.
✓ Branch 27 taken 2 times.
✓ Branch 28 taken 2 times.
✗ Branch 29 not taken.
✗ Branch 30 not taken.
✓ Branch 31 taken 2 times.
✓ Branch 32 taken 2 times.
✗ Branch 33 not taken.
✓ Branch 34 taken 2 times.
✗ Branch 35 not taken.
✗ Branch 36 not taken.
✓ Branch 37 taken 2 times.
✗ Branch 38 not taken.
✗ Branch 39 not taken.
✗ Branch 40 not taken.
✓ Branch 41 taken 2 times.
✗ Branch 42 not taken.
✗ Branch 43 not taken.
✗ Branch 44 not taken.
✓ Branch 45 taken 2 times.
✓ Branch 46 taken 2 times.
✗ Branch 47 not taken.
✗ Branch 48 not taken.
✓ Branch 49 taken 2 times.
✓ Branch 50 taken 2 times.
✗ Branch 51 not taken.
✓ Branch 52 taken 2 times.
✗ Branch 53 not taken.
✗ Branch 54 not taken.
✓ Branch 55 taken 2 times.
✗ Branch 56 not taken.
✗ Branch 57 not taken.
✗ Branch 58 not taken.
✓ Branch 59 taken 2 times.
✗ Branch 60 not taken.
✗ Branch 61 not taken.
|
136 | PetscValidLogicalCollectiveEnum(nep,refine,2); |
| 965 |
27/62✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✓ Branch 13 taken 2 times.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
✓ Branch 16 taken 2 times.
✗ Branch 17 not taken.
✗ Branch 18 not taken.
✓ Branch 19 taken 2 times.
✓ Branch 20 taken 2 times.
✗ Branch 21 not taken.
✗ Branch 22 not taken.
✓ Branch 23 taken 2 times.
✓ Branch 24 taken 2 times.
✗ Branch 25 not taken.
✗ Branch 26 not taken.
✓ Branch 27 taken 2 times.
✓ Branch 28 taken 2 times.
✗ Branch 29 not taken.
✗ Branch 30 not taken.
✓ Branch 31 taken 2 times.
✓ Branch 32 taken 2 times.
✗ Branch 33 not taken.
✓ Branch 34 taken 2 times.
✗ Branch 35 not taken.
✗ Branch 36 not taken.
✓ Branch 37 taken 2 times.
✗ Branch 38 not taken.
✗ Branch 39 not taken.
✗ Branch 40 not taken.
✓ Branch 41 taken 2 times.
✗ Branch 42 not taken.
✗ Branch 43 not taken.
✗ Branch 44 not taken.
✓ Branch 45 taken 2 times.
✓ Branch 46 taken 2 times.
✗ Branch 47 not taken.
✗ Branch 48 not taken.
✓ Branch 49 taken 2 times.
✓ Branch 50 taken 2 times.
✗ Branch 51 not taken.
✓ Branch 52 taken 2 times.
✗ Branch 53 not taken.
✗ Branch 54 not taken.
✓ Branch 55 taken 2 times.
✗ Branch 56 not taken.
✗ Branch 57 not taken.
✗ Branch 58 not taken.
✓ Branch 59 taken 2 times.
✗ Branch 60 not taken.
✗ Branch 61 not taken.
|
136 | PetscValidLogicalCollectiveInt(nep,npart,3); |
| 966 |
29/66✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
✗ Branch 16 not taken.
✓ Branch 17 taken 2 times.
✓ Branch 18 taken 2 times.
✗ Branch 19 not taken.
✗ Branch 20 not taken.
✓ Branch 21 taken 2 times.
✓ Branch 22 taken 2 times.
✗ Branch 23 not taken.
✗ Branch 24 not taken.
✓ Branch 25 taken 2 times.
✓ Branch 26 taken 2 times.
✗ Branch 27 not taken.
✗ Branch 28 not taken.
✓ Branch 29 taken 2 times.
✓ Branch 30 taken 2 times.
✗ Branch 31 not taken.
✗ Branch 32 not taken.
✓ Branch 33 taken 2 times.
✓ Branch 34 taken 2 times.
✗ Branch 35 not taken.
✓ Branch 36 taken 2 times.
✗ Branch 37 not taken.
✗ Branch 38 not taken.
✓ Branch 39 taken 2 times.
✗ Branch 40 not taken.
✗ Branch 41 not taken.
✗ Branch 42 not taken.
✓ Branch 43 taken 2 times.
✗ Branch 44 not taken.
✗ Branch 45 not taken.
✗ Branch 46 not taken.
✓ Branch 47 taken 2 times.
✓ Branch 48 taken 2 times.
✗ Branch 49 not taken.
✗ Branch 50 not taken.
✓ Branch 51 taken 2 times.
✓ Branch 52 taken 2 times.
✗ Branch 53 not taken.
✓ Branch 54 taken 2 times.
✗ Branch 55 not taken.
✗ Branch 56 not taken.
✓ Branch 57 taken 2 times.
✗ Branch 58 not taken.
✗ Branch 59 not taken.
✓ Branch 60 taken 2 times.
✗ Branch 61 not taken.
✗ Branch 62 not taken.
✓ Branch 63 taken 2 times.
✗ Branch 64 not taken.
✗ Branch 65 not taken.
|
136 | PetscValidLogicalCollectiveReal(nep,tol,4); |
| 967 |
27/62✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✓ Branch 13 taken 2 times.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
✓ Branch 16 taken 2 times.
✗ Branch 17 not taken.
✗ Branch 18 not taken.
✓ Branch 19 taken 2 times.
✓ Branch 20 taken 2 times.
✗ Branch 21 not taken.
✗ Branch 22 not taken.
✓ Branch 23 taken 2 times.
✓ Branch 24 taken 2 times.
✗ Branch 25 not taken.
✗ Branch 26 not taken.
✓ Branch 27 taken 2 times.
✓ Branch 28 taken 2 times.
✗ Branch 29 not taken.
✗ Branch 30 not taken.
✓ Branch 31 taken 2 times.
✓ Branch 32 taken 2 times.
✗ Branch 33 not taken.
✓ Branch 34 taken 2 times.
✗ Branch 35 not taken.
✗ Branch 36 not taken.
✓ Branch 37 taken 2 times.
✗ Branch 38 not taken.
✗ Branch 39 not taken.
✗ Branch 40 not taken.
✓ Branch 41 taken 2 times.
✗ Branch 42 not taken.
✗ Branch 43 not taken.
✗ Branch 44 not taken.
✓ Branch 45 taken 2 times.
✓ Branch 46 taken 2 times.
✗ Branch 47 not taken.
✗ Branch 48 not taken.
✓ Branch 49 taken 2 times.
✓ Branch 50 taken 2 times.
✗ Branch 51 not taken.
✓ Branch 52 taken 2 times.
✗ Branch 53 not taken.
✗ Branch 54 not taken.
✓ Branch 55 taken 2 times.
✗ Branch 56 not taken.
✗ Branch 57 not taken.
✗ Branch 58 not taken.
✓ Branch 59 taken 2 times.
✗ Branch 60 not taken.
✗ Branch 61 not taken.
|
136 | PetscValidLogicalCollectiveInt(nep,its,5); |
| 968 |
27/62✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✓ Branch 13 taken 2 times.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
✓ Branch 16 taken 2 times.
✗ Branch 17 not taken.
✗ Branch 18 not taken.
✓ Branch 19 taken 2 times.
✓ Branch 20 taken 2 times.
✗ Branch 21 not taken.
✗ Branch 22 not taken.
✓ Branch 23 taken 2 times.
✓ Branch 24 taken 2 times.
✗ Branch 25 not taken.
✗ Branch 26 not taken.
✓ Branch 27 taken 2 times.
✓ Branch 28 taken 2 times.
✗ Branch 29 not taken.
✗ Branch 30 not taken.
✓ Branch 31 taken 2 times.
✓ Branch 32 taken 2 times.
✗ Branch 33 not taken.
✓ Branch 34 taken 2 times.
✗ Branch 35 not taken.
✗ Branch 36 not taken.
✓ Branch 37 taken 2 times.
✗ Branch 38 not taken.
✗ Branch 39 not taken.
✗ Branch 40 not taken.
✓ Branch 41 taken 2 times.
✗ Branch 42 not taken.
✗ Branch 43 not taken.
✗ Branch 44 not taken.
✓ Branch 45 taken 2 times.
✓ Branch 46 taken 2 times.
✗ Branch 47 not taken.
✗ Branch 48 not taken.
✓ Branch 49 taken 2 times.
✓ Branch 50 taken 2 times.
✗ Branch 51 not taken.
✓ Branch 52 taken 2 times.
✗ Branch 53 not taken.
✗ Branch 54 not taken.
✓ Branch 55 taken 2 times.
✗ Branch 56 not taken.
✗ Branch 57 not taken.
✗ Branch 58 not taken.
✓ Branch 59 taken 2 times.
✗ Branch 60 not taken.
✗ Branch 61 not taken.
|
136 | PetscValidLogicalCollectiveEnum(nep,scheme,6); |
| 969 | 136 | nep->refine = refine; | |
| 970 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
136 | if (refine) { /* process parameters only if not REFINE_NONE */ |
| 971 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
136 | if (npart!=nep->npart) { |
| 972 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
80 | PetscCall(PetscSubcommDestroy(&nep->refinesubc)); |
| 973 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
80 | PetscCall(KSPDestroy(&nep->refineksp)); |
| 974 | } | ||
| 975 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
|
136 | if (npart == PETSC_DETERMINE) { |
| 976 | ✗ | nep->npart = 1; | |
| 977 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
136 | } else if (npart != PETSC_CURRENT) { |
| 978 |
14/28✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
✓ Branch 16 taken 2 times.
✗ Branch 17 not taken.
✓ Branch 18 taken 2 times.
✗ Branch 19 not taken.
✗ Branch 20 not taken.
✓ Branch 21 taken 2 times.
✗ Branch 22 not taken.
✗ Branch 23 not taken.
✗ Branch 24 not taken.
✓ Branch 25 taken 2 times.
✗ Branch 26 not taken.
✗ Branch 27 not taken.
|
136 | PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)nep),&size)); |
| 979 |
2/6✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 10 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
|
136 | PetscCheck(npart>0 && npart<=size,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of npart"); |
| 980 | 136 | nep->npart = npart; | |
| 981 | } | ||
| 982 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
136 | if (tol == (PetscReal)PETSC_DETERMINE) { |
| 983 | 120 | nep->rtol = PETSC_DETERMINE; | |
| 984 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
16 | } else if (tol != (PetscReal)PETSC_CURRENT) { |
| 985 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
16 | PetscCheck(tol>0.0,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of tol. Must be > 0"); |
| 986 | 16 | nep->rtol = tol; | |
| 987 | } | ||
| 988 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
136 | if (its==PETSC_DETERMINE) { |
| 989 | 110 | nep->rits = PETSC_DETERMINE; | |
| 990 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
26 | } else if (its != PETSC_CURRENT) { |
| 991 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
26 | PetscCheck(its>=0,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of its. Must be >= 0"); |
| 992 | 26 | nep->rits = its; | |
| 993 | } | ||
| 994 | 136 | nep->scheme = scheme; | |
| 995 | } | ||
| 996 | 136 | nep->state = NEP_STATE_INITIAL; | |
| 997 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
136 | PetscFunctionReturn(PETSC_SUCCESS); |
| 998 | } | ||
| 999 | |||
| 1000 | /*@ | ||
| 1001 | NEPGetRefine - Gets the refinement strategy used by the `NEP` object, and the | ||
| 1002 | associated parameters. | ||
| 1003 | |||
| 1004 | Not Collective | ||
| 1005 | |||
| 1006 | Input Parameter: | ||
| 1007 | . nep - the nonlinear eigensolver context | ||
| 1008 | |||
| 1009 | Output Parameters: | ||
| 1010 | + refine - refinement type | ||
| 1011 | . npart - number of partitions of the communicator | ||
| 1012 | . tol - the convergence tolerance | ||
| 1013 | . its - maximum number of refinement iterations | ||
| 1014 | - scheme - the scheme used for solving linear systems | ||
| 1015 | |||
| 1016 | Level: intermediate | ||
| 1017 | |||
| 1018 | Note: | ||
| 1019 | The user can specify `NULL` for any parameter that is not needed. | ||
| 1020 | |||
| 1021 | .seealso: [](ch:nep), `NEPSetRefine()` | ||
| 1022 | @*/ | ||
| 1023 | 16 | PetscErrorCode NEPGetRefine(NEP nep,NEPRefine *refine,PetscInt *npart,PetscReal *tol,PetscInt *its,NEPRefineScheme *scheme) | |
| 1024 | { | ||
| 1025 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
16 | PetscFunctionBegin; |
| 1026 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
16 | PetscValidHeaderSpecific(nep,NEP_CLASSID,1); |
| 1027 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
16 | if (refine) *refine = nep->refine; |
| 1028 |
2/2✓ Branch 0 taken 6 times.
✓ Branch 1 taken 10 times.
|
16 | if (npart) *npart = nep->npart; |
| 1029 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
16 | if (tol) *tol = nep->rtol; |
| 1030 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
16 | if (its) *its = nep->rits; |
| 1031 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
16 | if (scheme) *scheme = nep->scheme; |
| 1032 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
16 | PetscFunctionReturn(PETSC_SUCCESS); |
| 1033 | } | ||
| 1034 | |||
| 1035 | /*@ | ||
| 1036 | NEPSetOptionsPrefix - Sets the prefix used for searching for all | ||
| 1037 | `NEP` options in the database. | ||
| 1038 | |||
| 1039 | Logically Collective | ||
| 1040 | |||
| 1041 | Input Parameters: | ||
| 1042 | + nep - the nonlinear eigensolver context | ||
| 1043 | - prefix - the prefix string to prepend to all `NEP` option requests | ||
| 1044 | |||
| 1045 | Notes: | ||
| 1046 | A hyphen (-) must NOT be given at the beginning of the prefix name. | ||
| 1047 | The first character of all runtime options is AUTOMATICALLY the | ||
| 1048 | hyphen. | ||
| 1049 | |||
| 1050 | For example, to distinguish between the runtime options for two | ||
| 1051 | different `NEP` contexts, one could call | ||
| 1052 | .vb | ||
| 1053 | NEPSetOptionsPrefix(nep1,"neig1_") | ||
| 1054 | NEPSetOptionsPrefix(nep2,"neig2_") | ||
| 1055 | .ve | ||
| 1056 | |||
| 1057 | Level: advanced | ||
| 1058 | |||
| 1059 | .seealso: [](ch:nep), `NEPAppendOptionsPrefix()`, `NEPGetOptionsPrefix()` | ||
| 1060 | @*/ | ||
| 1061 | 24 | PetscErrorCode NEPSetOptionsPrefix(NEP nep,const char prefix[]) | |
| 1062 | { | ||
| 1063 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
24 | PetscFunctionBegin; |
| 1064 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
24 | PetscValidHeaderSpecific(nep,NEP_CLASSID,1); |
| 1065 |
6/8✓ Branch 0 taken 10 times.
✓ Branch 1 taken 2 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
24 | if (!nep->V) PetscCall(NEPGetBV(nep,&nep->V)); |
| 1066 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
24 | PetscCall(BVSetOptionsPrefix(nep->V,prefix)); |
| 1067 |
6/8✓ Branch 0 taken 10 times.
✓ Branch 1 taken 2 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
24 | if (!nep->ds) PetscCall(NEPGetDS(nep,&nep->ds)); |
| 1068 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
24 | PetscCall(DSSetOptionsPrefix(nep->ds,prefix)); |
| 1069 |
6/8✓ Branch 0 taken 10 times.
✓ Branch 1 taken 2 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
24 | if (!nep->rg) PetscCall(NEPGetRG(nep,&nep->rg)); |
| 1070 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
24 | PetscCall(RGSetOptionsPrefix(nep->rg,prefix)); |
| 1071 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
24 | PetscCall(PetscObjectSetOptionsPrefix((PetscObject)nep,prefix)); |
| 1072 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
2 | PetscFunctionReturn(PETSC_SUCCESS); |
| 1073 | } | ||
| 1074 | |||
| 1075 | /*@ | ||
| 1076 | NEPAppendOptionsPrefix - Appends to the prefix used for searching for all | ||
| 1077 | `NEP` options in the database. | ||
| 1078 | |||
| 1079 | Logically Collective | ||
| 1080 | |||
| 1081 | Input Parameters: | ||
| 1082 | + nep - the nonlinear eigensolver context | ||
| 1083 | - prefix - the prefix string to prepend to all `NEP` option requests | ||
| 1084 | |||
| 1085 | Notes: | ||
| 1086 | A hyphen (-) must NOT be given at the beginning of the prefix name. | ||
| 1087 | The first character of all runtime options is AUTOMATICALLY the hyphen. | ||
| 1088 | |||
| 1089 | Level: advanced | ||
| 1090 | |||
| 1091 | .seealso: [](ch:nep), `NEPSetOptionsPrefix()`, `NEPGetOptionsPrefix()` | ||
| 1092 | @*/ | ||
| 1093 | 10 | PetscErrorCode NEPAppendOptionsPrefix(NEP nep,const char prefix[]) | |
| 1094 | { | ||
| 1095 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
10 | PetscFunctionBegin; |
| 1096 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
10 | PetscValidHeaderSpecific(nep,NEP_CLASSID,1); |
| 1097 |
1/8✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
10 | if (!nep->V) PetscCall(NEPGetBV(nep,&nep->V)); |
| 1098 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(BVAppendOptionsPrefix(nep->V,prefix)); |
| 1099 |
1/8✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
10 | if (!nep->ds) PetscCall(NEPGetDS(nep,&nep->ds)); |
| 1100 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(DSAppendOptionsPrefix(nep->ds,prefix)); |
| 1101 |
1/8✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
10 | if (!nep->rg) PetscCall(NEPGetRG(nep,&nep->rg)); |
| 1102 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(RGAppendOptionsPrefix(nep->rg,prefix)); |
| 1103 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)nep,prefix)); |
| 1104 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
2 | PetscFunctionReturn(PETSC_SUCCESS); |
| 1105 | } | ||
| 1106 | |||
| 1107 | /*@ | ||
| 1108 | NEPGetOptionsPrefix - Gets the prefix used for searching for all | ||
| 1109 | `NEP` options in the database. | ||
| 1110 | |||
| 1111 | Not Collective | ||
| 1112 | |||
| 1113 | Input Parameter: | ||
| 1114 | . nep - the nonlinear eigensolver context | ||
| 1115 | |||
| 1116 | Output Parameter: | ||
| 1117 | . prefix - pointer to the prefix string used is returned | ||
| 1118 | |||
| 1119 | Level: advanced | ||
| 1120 | |||
| 1121 | .seealso: [](ch:nep), `NEPSetOptionsPrefix()`, `NEPAppendOptionsPrefix()` | ||
| 1122 | @*/ | ||
| 1123 | 31 | PetscErrorCode NEPGetOptionsPrefix(NEP nep,const char *prefix[]) | |
| 1124 | { | ||
| 1125 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
31 | PetscFunctionBegin; |
| 1126 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
31 | PetscValidHeaderSpecific(nep,NEP_CLASSID,1); |
| 1127 |
2/8✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
31 | PetscAssertPointer(prefix,2); |
| 1128 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
31 | PetscCall(PetscObjectGetOptionsPrefix((PetscObject)nep,prefix)); |
| 1129 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
2 | PetscFunctionReturn(PETSC_SUCCESS); |
| 1130 | } | ||
| 1131 |