Actual source code: nepsetup.c
1: /*
2: NEP routines related to problem setup.
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2013, 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*/
25: #include <slepc-private/ipimpl.h>
29: /*@
30: NEPSetUp - Sets up all the internal data structures necessary for the
31: execution of the NEP solver.
33: Collective on NEP
35: Input Parameter:
36: . nep - solver context
38: Notes:
39: This function need not be called explicitly in most cases, since NEPSolve()
40: calls it. It can be useful when one wants to measure the set-up time
41: separately from the solve time.
43: Level: advanced
45: .seealso: NEPCreate(), NEPSolve(), NEPDestroy()
46: @*/
47: PetscErrorCode NEPSetUp(NEP nep)
48: {
50: Mat T;
54: if (nep->setupcalled) return(0);
55: PetscLogEventBegin(NEP_SetUp,nep,0,0,0);
57: /* reset the convergence flag from the previous solves */
58: nep->reason = NEP_CONVERGED_ITERATING;
60: /* Set default solver type (NEPSetFromOptions was not called) */
61: if (!((PetscObject)nep)->type_name) {
62: NEPSetType(nep,NEPRII);
63: }
64: if (!nep->ip) { NEPGetIP(nep,&nep->ip); }
65: if (!((PetscObject)nep->ip)->type_name) {
66: IPSetType_Default(nep->ip);
67: }
68: if (!nep->ds) { NEPGetDS(nep,&nep->ds); }
69: DSReset(nep->ds);
70: if (!((PetscObject)nep->rand)->type_name) {
71: PetscRandomSetFromOptions(nep->rand);
72: }
73: if (!nep->ksp) {
74: NEPGetKSP(nep,&nep->ksp);
75: }
77: /* by default, compute eigenvalues close to target */
78: /* nep->target should contain the initial guess for the eigenvalue */
79: if (!nep->which) nep->which = NEP_TARGET_MAGNITUDE;
81: /* set problem dimensions */
82: VecDestroy(&nep->t);
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(nep,nep->function);
87: PetscLogObjectParent(nep,nep->jacobian);
88: SlepcMatGetVecsTemplate(nep->A[0],&nep->t,NULL);
89: } else {
90: NEPGetFunction(nep,&T,NULL,NULL,NULL);
91: SlepcMatGetVecsTemplate(T,&nep->t,NULL);
92: }
93: PetscLogObjectParent(nep,nep->t);
94: VecGetSize(nep->t,&nep->n);
95: VecGetLocalSize(nep->t,&nep->nloc);
97: /* call specific solver setup */
98: (*nep->ops->setup)(nep);
100: /* set tolerances if not yet set */
101: if (nep->abstol==PETSC_DEFAULT) nep->abstol = 1e-50;
102: if (nep->rtol==PETSC_DEFAULT) nep->rtol = 100*SLEPC_DEFAULT_TOL;
103: if (nep->stol==PETSC_DEFAULT) nep->stol = SLEPC_DEFAULT_TOL;
104: nep->ktol = 0.1;
105: nep->nfuncs = 0;
106: nep->linits = 0;
108: /* set eigenvalue comparison */
109: switch (nep->which) {
110: case NEP_LARGEST_MAGNITUDE:
111: nep->comparison = SlepcCompareLargestMagnitude;
112: nep->comparisonctx = NULL;
113: break;
114: case NEP_SMALLEST_MAGNITUDE:
115: nep->comparison = SlepcCompareSmallestMagnitude;
116: nep->comparisonctx = NULL;
117: break;
118: case NEP_LARGEST_REAL:
119: nep->comparison = SlepcCompareLargestReal;
120: nep->comparisonctx = NULL;
121: break;
122: case NEP_SMALLEST_REAL:
123: nep->comparison = SlepcCompareSmallestReal;
124: nep->comparisonctx = NULL;
125: break;
126: case NEP_LARGEST_IMAGINARY:
127: nep->comparison = SlepcCompareLargestImaginary;
128: nep->comparisonctx = NULL;
129: break;
130: case NEP_SMALLEST_IMAGINARY:
131: nep->comparison = SlepcCompareSmallestImaginary;
132: nep->comparisonctx = NULL;
133: break;
134: case NEP_TARGET_MAGNITUDE:
135: nep->comparison = SlepcCompareTargetMagnitude;
136: nep->comparisonctx = &nep->target;
137: break;
138: case NEP_TARGET_REAL:
139: nep->comparison = SlepcCompareTargetReal;
140: nep->comparisonctx = &nep->target;
141: break;
142: case NEP_TARGET_IMAGINARY:
143: nep->comparison = SlepcCompareTargetImaginary;
144: nep->comparisonctx = &nep->target;
145: break;
146: }
148: if (nep->ncv > 2*nep->n) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"ncv must be twice the problem size at most");
149: if (nep->nev > nep->ncv) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"nev bigger than ncv");
151: /* process initial vectors */
152: if (nep->nini<0) {
153: nep->nini = -nep->nini;
154: if (nep->nini>nep->ncv) SETERRQ(PetscObjectComm((PetscObject)nep),1,"The number of initial vectors is larger than ncv");
155: IPOrthonormalizeBasis_Private(nep->ip,&nep->nini,&nep->IS,nep->V);
156: }
157: PetscLogEventEnd(NEP_SetUp,nep,0,0,0);
158: nep->setupcalled = 1;
159: return(0);
160: }
164: /*@
165: NEPSetInitialSpace - Specify a basis of vectors that constitute the initial
166: space, that is, the subspace from which the solver starts to iterate.
168: Collective on NEP and Vec
170: Input Parameter:
171: + nep - the nonlinear eigensolver context
172: . n - number of vectors
173: - is - set of basis vectors of the initial space
175: Notes:
176: Some solvers start to iterate on a single vector (initial vector). In that case,
177: the other vectors are ignored.
179: These vectors do not persist from one NEPSolve() call to the other, so the
180: initial space should be set every time.
182: The vectors do not need to be mutually orthonormal, since they are explicitly
183: orthonormalized internally.
185: Common usage of this function is when the user can provide a rough approximation
186: of the wanted eigenspace. Then, convergence may be faster.
188: Level: intermediate
189: @*/
190: PetscErrorCode NEPSetInitialSpace(NEP nep,PetscInt n,Vec *is)
191: {
197: if (n<0) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Argument n cannot be negative");
198: SlepcBasisReference_Private(n,is,&nep->nini,&nep->IS);
199: if (n>0) nep->setupcalled = 0;
200: return(0);
201: }
205: /*
206: NEPAllocateSolution - Allocate memory storage for common variables such
207: as eigenvalues and eigenvectors. All vectors in V (and W) share a
208: contiguous chunk of memory.
209: */
210: PetscErrorCode NEPAllocateSolution(NEP nep)
211: {
213: PetscInt newc,cnt;
216: if (nep->allocated_ncv != nep->ncv) {
217: newc = PetscMax(0,nep->ncv-nep->allocated_ncv);
218: NEPFreeSolution(nep);
219: cnt = 0;
220: PetscMalloc(nep->ncv*sizeof(PetscScalar),&nep->eig);
221: cnt += 2*newc*sizeof(PetscScalar);
222: PetscMalloc(nep->ncv*sizeof(PetscReal),&nep->errest);
223: PetscMalloc(nep->ncv*sizeof(PetscInt),&nep->perm);
224: cnt += 2*newc*sizeof(PetscReal);
225: PetscLogObjectMemory(nep,cnt);
226: VecDuplicateVecs(nep->t,nep->ncv,&nep->V);
227: PetscLogObjectParents(nep,nep->ncv,nep->V);
228: nep->allocated_ncv = nep->ncv;
229: }
230: return(0);
231: }
235: /*
236: NEPFreeSolution - Free memory storage. This routine is related to
237: NEPAllocateSolution().
238: */
239: PetscErrorCode NEPFreeSolution(NEP nep)
240: {
244: if (nep->allocated_ncv > 0) {
245: PetscFree(nep->eig);
246: PetscFree(nep->errest);
247: PetscFree(nep->perm);
248: VecDestroyVecs(nep->allocated_ncv,&nep->V);
249: nep->allocated_ncv = 0;
250: }
251: return(0);
252: }