Actual source code: svdsetup.c

  1: /*
  2:      SVD routines for setting up the solver.

  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/svdimpl.h>      /*I "slepcsvd.h" I*/
 25: #include <slepc-private/ipimpl.h>

 29: /*@
 30:    SVDSetOperator - Set the matrix associated with the singular value problem.

 32:    Collective on SVD and Mat

 34:    Input Parameters:
 35: +  svd - the singular value solver context
 36: -  A  - the matrix associated with the singular value problem

 38:    Level: beginner

 40: .seealso: SVDSolve(), SVDGetOperator()
 41: @*/
 42: PetscErrorCode SVDSetOperator(SVD svd,Mat mat)
 43: {

 50:   if (svd->setupcalled) { SVDReset(svd); }
 51:   PetscObjectReference((PetscObject)mat);
 52:   MatDestroy(&svd->OP);
 53:   svd->OP = mat;
 54:   return(0);
 55: }

 59: /*@
 60:    SVDGetOperator - Get the matrix associated with the singular value problem.

 62:    Not collective, though parallel Mats are returned if the SVD is parallel

 64:    Input Parameter:
 65: .  svd - the singular value solver context

 67:    Output Parameters:
 68: .  A    - the matrix associated with the singular value problem

 70:    Level: advanced

 72: .seealso: SVDSolve(), SVDSetOperator()
 73: @*/
 74: PetscErrorCode SVDGetOperator(SVD svd,Mat *A)
 75: {
 79:   *A = svd->OP;
 80:   return(0);
 81: }

 85: /*@
 86:    SVDSetUp - Sets up all the internal data structures necessary for the
 87:    execution of the singular value solver.

 89:    Collective on SVD

 91:    Input Parameter:
 92: .  svd   - singular value solver context

 94:    Level: advanced

 96:    Notes:
 97:    This function need not be called explicitly in most cases, since SVDSolve()
 98:    calls it. It can be useful when one wants to measure the set-up time
 99:    separately from the solve time.

101: .seealso: SVDCreate(), SVDSolve(), SVDDestroy()
102: @*/
103: PetscErrorCode SVDSetUp(SVD svd)
104: {
106:   PetscBool      flg;
107:   PetscInt       M,N,k;
108:   Vec            *T;

112:   if (svd->setupcalled) return(0);
113:   PetscLogEventBegin(SVD_SetUp,svd,0,0,0);

115:   /* reset the convergence flag from the previous solves */
116:   svd->reason = SVD_CONVERGED_ITERATING;

118:   /* Set default solver type (SVDSetFromOptions was not called) */
119:   if (!((PetscObject)svd)->type_name) {
120:     SVDSetType(svd,SVDCROSS);
121:   }
122:   if (!svd->ip) { SVDGetIP(svd,&svd->ip); }
123:   if (!((PetscObject)svd->ip)->type_name) {
124:     IPSetType_Default(svd->ip);
125:   }
126:   if (!svd->ds) { SVDGetDS(svd,&svd->ds); }
127:   DSReset(svd->ds);
128:   if (!((PetscObject)svd->rand)->type_name) {
129:     PetscRandomSetFromOptions(svd->rand);
130:   }

132:   /* check matrix */
133:   if (!svd->OP) SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_WRONGSTATE,"SVDSetOperator must be called first");

135:   /* determine how to build the transpose */
136:   if (svd->transmode == PETSC_DECIDE) {
137:     MatHasOperation(svd->OP,MATOP_TRANSPOSE,&flg);
138:     if (flg) svd->transmode = SVD_TRANSPOSE_EXPLICIT;
139:     else svd->transmode = SVD_TRANSPOSE_IMPLICIT;
140:   }

142:   /* build transpose matrix */
143:   MatDestroy(&svd->A);
144:   MatDestroy(&svd->AT);
145:   MatGetSize(svd->OP,&M,&N);
146:   PetscObjectReference((PetscObject)svd->OP);
147:   switch (svd->transmode) {
148:     case SVD_TRANSPOSE_EXPLICIT:
149:       if (M>=N) {
150:         svd->A = svd->OP;
151:         MatTranspose(svd->OP,MAT_INITIAL_MATRIX,&svd->AT);
152:         MatConjugate(svd->AT);
153:       } else {
154:         MatTranspose(svd->OP,MAT_INITIAL_MATRIX,&svd->A);
155:         MatConjugate(svd->A);
156:         svd->AT = svd->OP;
157:       }
158:       break;
159:     case SVD_TRANSPOSE_IMPLICIT:
160:       if (M>=N) {
161:         svd->A = svd->OP;
162:         svd->AT = NULL;
163:       } else {
164:         svd->A = NULL;
165:         svd->AT = svd->OP;
166:       }
167:       break;
168:     default:
169:       SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Invalid transpose mode");
170:   }

172:   VecDestroy(&svd->tr);
173:   VecDestroy(&svd->tl);
174:   if (svd->A) {
175:     SlepcMatGetVecsTemplate(svd->A,&svd->tr,&svd->tl);
176:   } else {
177:     SlepcMatGetVecsTemplate(svd->AT,&svd->tl,&svd->tr);
178:   }
179:   PetscLogObjectParent(svd,svd->tl);
180:   PetscLogObjectParent(svd,svd->tr);

182:   /* swap initial vectors if necessary */
183:   if (M<N) {
184:     T=svd->ISL; svd->ISL=svd->IS; svd->IS=T;
185:     k=svd->ninil; svd->ninil=svd->nini; svd->nini=k;
186:   }

188:   /* call specific solver setup */
189:   (*svd->ops->setup)(svd);

191:   /* set tolerance if not yet set */
192:   if (svd->tol==PETSC_DEFAULT) svd->tol = SLEPC_DEFAULT_TOL;

194:   if (svd->ncv > M || svd->ncv > N) SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"ncv bigger than matrix dimensions");
195:   if (svd->nsv > svd->ncv) SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"nsv bigger than ncv");

197:   if (svd->ncv != svd->n) {
198:     /* free memory for previous solution  */
199:     if (svd->n) {
200:       PetscFree(svd->sigma);
201:       PetscFree(svd->perm);
202:       PetscFree(svd->errest);
203:       VecDestroyVecs(svd->n,&svd->V);
204:     }
205:     /* allocate memory for next solution */
206:     PetscMalloc(svd->ncv*sizeof(PetscReal),&svd->sigma);
207:     PetscMalloc(svd->ncv*sizeof(PetscInt),&svd->perm);
208:     PetscMalloc(svd->ncv*sizeof(PetscReal),&svd->errest);
209:     PetscLogObjectMemory(svd,PetscMax(0,svd->ncv-svd->n)*(2*sizeof(PetscReal)+sizeof(PetscInt)));
210:     VecDuplicateVecs(svd->tr,svd->ncv,&svd->V);
211:     PetscLogObjectParents(svd,svd->ncv,svd->V);
212:     svd->n = svd->ncv;
213:   }

215:   /* process initial vectors */
216:   if (svd->nini<0) {
217:     svd->nini = -svd->nini;
218:     if (svd->nini>svd->ncv) SETERRQ(PetscObjectComm((PetscObject)svd),1,"The number of initial vectors is larger than ncv");
219:     IPOrthonormalizeBasis_Private(svd->ip,&svd->nini,&svd->IS,svd->V);
220:   }
221:   if (svd->ninil<0 && svd->U) { /* skip this if the solver is not using a left basis */
222:     svd->ninil = -svd->ninil;
223:     if (svd->ninil>svd->ncv) SETERRQ(PetscObjectComm((PetscObject)svd),1,"The number of left initial vectors is larger than ncv");
224:     IPOrthonormalizeBasis_Private(svd->ip,&svd->ninil,&svd->ISL,svd->U);
225:   }

227:   PetscLogEventEnd(SVD_SetUp,svd,0,0,0);
228:   svd->setupcalled = 1;
229:   return(0);
230: }

234: /*@
235:    SVDSetInitialSpace - Specify a basis of vectors that constitute the initial
236:    (right) space, that is, a rough approximation to the right singular subspace
237:    from which the solver starts to iterate.

239:    Collective on SVD and Vec

241:    Input Parameter:
242: +  svd   - the singular value solver context
243: .  n     - number of vectors
244: -  is    - set of basis vectors of the initial space

246:    Notes:
247:    Some solvers start to iterate on a single vector (initial vector). In that case,
248:    the other vectors are ignored.

250:    These vectors do not persist from one SVDSolve() call to the other, so the
251:    initial space should be set every time.

253:    The vectors do not need to be mutually orthonormal, since they are explicitly
254:    orthonormalized internally.

256:    Common usage of this function is when the user can provide a rough approximation
257:    of the wanted singular space. Then, convergence may be faster.

259:    Level: intermediate
260: @*/
261: PetscErrorCode SVDSetInitialSpace(SVD svd,PetscInt n,Vec *is)
262: {

268:   if (n<0) SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Argument n cannot be negative");
269:   SlepcBasisReference_Private(n,is,&svd->nini,&svd->IS);
270:   if (n>0) svd->setupcalled = 0;
271:   return(0);
272: }

276: /*@
277:    SVDSetInitialSpaceLeft - Specify a basis of vectors that constitute the initial
278:    left space, that is, a rough approximation to the left singular subspace
279:    from which the solver starts to iterate.

281:    Collective on SVD and Vec

283:    Input Parameter:
284: +  svd   - the singular value solver context
285: .  n     - number of vectors
286: -  is    - set of basis vectors of the initial space

288:    Notes:
289:    Some solvers start to iterate on a single vector (initial vector). In that case,
290:    the other vectors are ignored.

292:    These vectors do not persist from one SVDSolve() call to the other, so the
293:    initial space should be set every time.

295:    The vectors do not need to be mutually orthonormal, since they are explicitly
296:    orthonormalized internally.

298:    Common usage of this function is when the user can provide a rough approximation
299:    of the wanted singular space. Then, convergence may be faster.

301:    Level: intermediate
302: @*/
303: PetscErrorCode SVDSetInitialSpaceLeft(SVD svd,PetscInt n,Vec *is)
304: {

310:   if (n<0) SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Argument n cannot be negative");
311:   SlepcBasisReference_Private(n,is,&svd->ninil,&svd->ISL);
312:   if (n>0) svd->setupcalled = 0;
313:   return(0);
314: }