1: /*
2: NEP routines related to problem setup.
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2015, Universitat Politecnica de Valencia, Spain
8: This file is part of SLEPc.
10: SLEPc is free software: you can redistribute it and/or modify it under the
11: terms of version 3 of the GNU Lesser General Public License as published by
12: the Free Software Foundation.
14: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
15: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
16: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
17: more details.
19: You should have received a copy of the GNU Lesser General Public License
20: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
21: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
22: */
24: #include <slepc/private/nepimpl.h> /*I "slepcnep.h" I*/
28: /*@
29: NEPSetUp - Sets up all the internal data structures necessary for the
30: execution of the NEP solver.
32: Collective on NEP 34: Input Parameter:
35: . nep - solver context
37: Notes:
38: This function need not be called explicitly in most cases, since NEPSolve()
39: calls it. It can be useful when one wants to measure the set-up time
40: separately from the solve time.
42: Level: developer
44: .seealso: NEPCreate(), NEPSolve(), NEPDestroy()
45: @*/
46: PetscErrorCode NEPSetUp(NEP nep) 47: {
49: PetscInt k;
50: SlepcSC sc;
51: Mat T;
55: if (nep->state) return(0);
56: PetscLogEventBegin(NEP_SetUp,nep,0,0,0);
58: /* reset the convergence flag from the previous solves */
59: nep->reason = NEP_CONVERGED_ITERATING;
61: /* set default solver type (NEPSetFromOptions was not called) */
62: if (!((PetscObject)nep)->type_name) {
63: NEPSetType(nep,NEPRII);
64: }
65: if (!nep->ds) { NEPGetDS(nep,&nep->ds); }
66: DSReset(nep->ds);
67: if (!nep->rg) { NEPGetRG(nep,&nep->rg); }
68: if (!((PetscObject)nep->rg)->type_name) {
69: RGSetType(nep->rg,RGINTERVAL);
70: }
71: if (!((PetscObject)nep->rand)->type_name) {
72: PetscRandomSetFromOptions(nep->rand);
73: }
74: if (!nep->ksp) {
75: NEPGetKSP(nep,&nep->ksp);
76: }
78: /* by default, compute eigenvalues close to target */
79: /* nep->target should contain the initial guess for the eigenvalue */
80: if (!nep->which) nep->which = NEP_TARGET_MAGNITUDE;
82: /* set problem dimensions */
83: if (nep->split) {
84: MatDuplicate(nep->A[0],MAT_DO_NOT_COPY_VALUES,&nep->function);
85: MatDuplicate(nep->A[0],MAT_DO_NOT_COPY_VALUES,&nep->jacobian);
86: PetscLogObjectParent((PetscObject)nep,(PetscObject)nep->function);
87: PetscLogObjectParent((PetscObject)nep,(PetscObject)nep->jacobian);
88: MatGetSize(nep->A[0],&nep->n,NULL);
89: MatGetLocalSize(nep->A[0],&nep->nloc,NULL);
90: } else {
91: NEPGetFunction(nep,&T,NULL,NULL,NULL);
92: MatGetSize(T,&nep->n,NULL);
93: MatGetLocalSize(T,&nep->nloc,NULL);
94: }
96: /* call specific solver setup */
97: (*nep->ops->setup)(nep);
99: /* set tolerances if not yet set */
100: if (nep->abstol==PETSC_DEFAULT) nep->abstol = 1e-50;
101: if (nep->rtol==PETSC_DEFAULT) nep->rtol = 100*SLEPC_DEFAULT_TOL;
102: if (nep->stol==PETSC_DEFAULT) nep->stol = SLEPC_DEFAULT_TOL;
103: nep->ktol = 0.1;
104: nep->nfuncs = 0;
105: if (nep->refine) {
106: if (nep->reftol==PETSC_DEFAULT) nep->reftol = SLEPC_DEFAULT_TOL;
107: if (nep->rits==PETSC_DEFAULT) nep->rits = (nep->refine==NEP_REFINE_SIMPLE)? 10: 1;
108: }
110: /* fill sorting criterion context */
111: switch (nep->which) {
112: case NEP_LARGEST_MAGNITUDE:
113: nep->sc->comparison = SlepcCompareLargestMagnitude;
114: nep->sc->comparisonctx = NULL;
115: break;
116: case NEP_SMALLEST_MAGNITUDE:
117: nep->sc->comparison = SlepcCompareSmallestMagnitude;
118: nep->sc->comparisonctx = NULL;
119: break;
120: case NEP_LARGEST_REAL:
121: nep->sc->comparison = SlepcCompareLargestReal;
122: nep->sc->comparisonctx = NULL;
123: break;
124: case NEP_SMALLEST_REAL:
125: nep->sc->comparison = SlepcCompareSmallestReal;
126: nep->sc->comparisonctx = NULL;
127: break;
128: case NEP_LARGEST_IMAGINARY:
129: nep->sc->comparison = SlepcCompareLargestImaginary;
130: nep->sc->comparisonctx = NULL;
131: break;
132: case NEP_SMALLEST_IMAGINARY:
133: nep->sc->comparison = SlepcCompareSmallestImaginary;
134: nep->sc->comparisonctx = NULL;
135: break;
136: case NEP_TARGET_MAGNITUDE:
137: nep->sc->comparison = SlepcCompareTargetMagnitude;
138: nep->sc->comparisonctx = &nep->target;
139: break;
140: case NEP_TARGET_REAL:
141: nep->sc->comparison = SlepcCompareTargetReal;
142: nep->sc->comparisonctx = &nep->target;
143: break;
144: case NEP_TARGET_IMAGINARY:
145: nep->sc->comparison = SlepcCompareTargetImaginary;
146: nep->sc->comparisonctx = &nep->target;
147: break;
148: }
150: nep->sc->map = NULL;
151: nep->sc->mapobj = NULL;
153: /* fill sorting criterion for DS */
154: DSGetSlepcSC(nep->ds,&sc);
155: sc->comparison = nep->sc->comparison;
156: sc->comparisonctx = nep->sc->comparisonctx;
157: sc->map = NULL;
158: sc->mapobj = NULL;
160: if (nep->ncv > nep->n) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"ncv must be the problem size at most");
161: if (nep->nev > nep->ncv) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"nev bigger than ncv");
163: /* process initial vectors */
164: if (nep->nini<0) {
165: k = -nep->nini;
166: if (k>nep->ncv) SETERRQ(PetscObjectComm((PetscObject)nep),1,"The number of initial vectors is larger than ncv");
167: BVInsertVecs(nep->V,0,&k,nep->IS,PETSC_TRUE);
168: SlepcBasisDestroy_Private(&nep->nini,&nep->IS);
169: nep->nini = k;
170: }
171: PetscLogEventEnd(NEP_SetUp,nep,0,0,0);
172: nep->state = NEP_STATE_SETUP;
173: return(0);
174: }
178: /*@
179: NEPSetInitialSpace - Specify a basis of vectors that constitute the initial
180: space, that is, the subspace from which the solver starts to iterate.
182: Collective on NEP and Vec
184: Input Parameter:
185: + nep - the nonlinear eigensolver context
186: . n - number of vectors
187: - is - set of basis vectors of the initial space
189: Notes:
190: Some solvers start to iterate on a single vector (initial vector). In that case,
191: the other vectors are ignored.
193: These vectors do not persist from one NEPSolve() call to the other, so the
194: initial space should be set every time.
196: The vectors do not need to be mutually orthonormal, since they are explicitly
197: orthonormalized internally.
199: Common usage of this function is when the user can provide a rough approximation
200: of the wanted eigenspace. Then, convergence may be faster.
202: Level: intermediate
203: @*/
204: PetscErrorCode NEPSetInitialSpace(NEP nep,PetscInt n,Vec *is)205: {
211: if (n<0) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Argument n cannot be negative");
212: SlepcBasisReference_Private(n,is,&nep->nini,&nep->IS);
213: if (n>0) nep->state = NEP_STATE_INITIAL;
214: return(0);
215: }
219: /*@
220: NEPAllocateSolution - Allocate memory storage for common variables such
221: as eigenvalues and eigenvectors.
223: Collective on NEP225: Input Parameters:
226: + nep - eigensolver context
227: - extra - number of additional positions, used for methods that require a
228: working basis slightly larger than ncv
230: Developers Note:
231: This is PETSC_EXTERN because it may be required by user plugin NEP232: implementations.
234: Level: developer
235: @*/
236: PetscErrorCode NEPAllocateSolution(NEP nep,PetscInt extra)237: {
239: PetscInt oldsize,newc,requested;
240: PetscLogDouble cnt;
241: Mat T;
242: Vec t;
245: requested = nep->ncv + extra;
247: /* oldsize is zero if this is the first time setup is called */
248: BVGetSizes(nep->V,NULL,NULL,&oldsize);
249: newc = PetscMax(0,requested-oldsize);
251: /* allocate space for eigenvalues and friends */
252: if (requested != oldsize || !nep->eigr) {
253: if (oldsize) {
254: PetscFree4(nep->eigr,nep->eigi,nep->errest,nep->perm);
255: }
256: PetscMalloc4(requested,&nep->eigr,requested,&nep->eigi,requested,&nep->errest,requested,&nep->perm);
257: cnt = newc*sizeof(PetscScalar) + newc*sizeof(PetscReal) + newc*sizeof(PetscInt);
258: PetscLogObjectMemory((PetscObject)nep,cnt);
259: }
261: /* allocate V */
262: if (!nep->V) { NEPGetBV(nep,&nep->V); }
263: if (!oldsize) {
264: if (!((PetscObject)(nep->V))->type_name) {
265: BVSetType(nep->V,BVSVEC);
266: }
267: if (nep->split) T = nep->A[0];
268: else {
269: NEPGetFunction(nep,&T,NULL,NULL,NULL);
270: }
271: MatCreateVecs(T,&t,NULL);
272: BVSetSizesFromVec(nep->V,t,requested);
273: VecDestroy(&t);
274: } else {
275: BVResize(nep->V,requested,PETSC_FALSE);
276: }
277: return(0);
278: }