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-2010, Universidad Politecnica de Valencia, Spain

  8:    This file is part of SLEPc.
  9:       
 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 private/svdimpl.h

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

 31:    Collective on SVD and Mat

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

 37:    Level: beginner

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

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

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

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

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

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

 71:    Level: advanced

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

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

 90:    Collective on SVD

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

 95:    Level: advanced

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

102: .seealso: SVDCreate(), SVDSolve(), SVDDestroy()
103: @*/
104: PetscErrorCode SVDSetUp(SVD svd)
105: {
107:   PetscTruth     flg,lindep;
108:   PetscInt       i,k,M,N,nloc;
109:   PetscScalar    *pV;
110:   PetscReal      norm;
111: 
114:   if (svd->setupcalled) return(0);
115:   PetscLogEventBegin(SVD_SetUp,svd,0,0,0);

117:   /* Set default solver type */
118:   if (!((PetscObject)svd)->type_name) {
119:     SVDSetType(svd,SVDCROSS);
120:   }

122:   /* check matrix */
123:   if (!svd->OP)
124:     SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "SVDSetOperator must be called first");
125: 
126:   /* determine how to build the transpose */
127:   if (svd->transmode == PETSC_DECIDE) {
128:     MatHasOperation(svd->OP,MATOP_TRANSPOSE,&flg);
129:     if (flg) svd->transmode = SVD_TRANSPOSE_EXPLICIT;
130:     else svd->transmode = SVD_TRANSPOSE_IMPLICIT;
131:   }
132: 
133:   /* build transpose matrix */
134:   if (svd->A) { MatDestroy(svd->A); }
135:   if (svd->AT) { MatDestroy(svd->AT); }
136:   MatGetSize(svd->OP,&M,&N);
137:   PetscObjectReference((PetscObject)svd->OP);
138:   switch (svd->transmode) {
139:     case SVD_TRANSPOSE_EXPLICIT:
140:       MatHasOperation(svd->OP,MATOP_TRANSPOSE,&flg);
141:       if (!flg) SETERRQ(1,"Matrix has not defined the MatTranpose operation");
142:       if (M>=N) {
143:         svd->A = svd->OP;
144:         MatTranspose(svd->OP, MAT_INITIAL_MATRIX,&svd->AT);
145:       } else {
146:         MatTranspose(svd->OP, MAT_INITIAL_MATRIX,&svd->A);
147:         svd->AT = svd->OP;
148:       }
149:       break;
150:     case SVD_TRANSPOSE_IMPLICIT:
151:       MatHasOperation(svd->OP,MATOP_MULT_TRANSPOSE,&flg);
152:       if (!flg) SETERRQ(1,"Matrix has not defined the MatMultTranpose operation");
153:       if (M>=N) {
154:         svd->A = svd->OP;
155:         svd->AT = PETSC_NULL;
156:       } else {
157:         svd->A = PETSC_NULL;
158:         svd->AT = svd->OP;
159:       }
160:       break;
161:     default:
162:       SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Invalid transpose mode");
163:   }

165:   /* initialize the random number generator */
166:   PetscRandomCreate(((PetscObject)svd)->comm,&svd->rand);
167:   PetscRandomSetFromOptions(svd->rand);

169:   /* call specific solver setup */
170:   (*svd->ops->setup)(svd);

172:   if (svd->ncv > M || svd->ncv > N)
173:     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"ncv bigger than matrix dimensions");
174:   if (svd->nsv > svd->ncv)
175:     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"nsv bigger than ncv");

177:   if (svd->ncv != svd->n) {
178:     /* free memory for previous solution  */
179:     if (svd->n) {
180:       PetscFree(svd->sigma);
181:       PetscFree(svd->perm);
182:       PetscFree(svd->errest);
183:       VecGetArray(svd->V[0],&pV);
184:       for (i=0;i<svd->n;i++) {
185:         VecDestroy(svd->V[i]);
186:       }
187:       PetscFree(pV);
188:       PetscFree(svd->V);
189:     }
190:     /* allocate memory for next solution */
191:     PetscMalloc(svd->ncv*sizeof(PetscReal),&svd->sigma);
192:     PetscMalloc(svd->ncv*sizeof(PetscInt),&svd->perm);
193:     PetscMalloc(svd->ncv*sizeof(PetscReal),&svd->errest);
194:     PetscMalloc(svd->ncv*sizeof(Vec),&svd->V);
195:     if (svd->A) {
196:       MatGetLocalSize(svd->A,PETSC_NULL,&nloc);
197:     } else {
198:       MatGetLocalSize(svd->AT,&nloc,PETSC_NULL);
199:     }
200:     PetscMalloc(svd->ncv*nloc*sizeof(PetscScalar),&pV);
201:     for (i=0;i<svd->ncv;i++) {
202:       VecCreateMPIWithArray(((PetscObject)svd)->comm,nloc,PETSC_DECIDE,pV+i*nloc,&svd->V[i]);
203:     }
204:     svd->n = svd->ncv;
205:   }

207:   /* process initial vectors */
208:   if (svd->nini<0) {
209:     svd->nini = -svd->nini;
210:     if (svd->nini>svd->ncv) SETERRQ(1,"The number of initial vectors is larger than ncv")
211:     k = 0;
212:     for (i=0;i<svd->nini;i++) {
213:       VecCopy(svd->IS[i],svd->V[k]);
214:       VecDestroy(svd->IS[i]);
215:       IPOrthogonalize(svd->ip,0,PETSC_NULL,k,PETSC_NULL,svd->V,svd->V[k],PETSC_NULL,&norm,&lindep);
216:       if (norm==0.0 || lindep) PetscInfo(svd,"Linearly dependent initial vector found, removing...\n");
217:       else {
218:         VecScale(svd->V[k],1.0/norm);
219:         k++;
220:       }
221:     }
222:     svd->nini = k;
223:     PetscFree(svd->IS);
224:   }

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

233: /*@
234:    SVDSetInitialSpace - Specify a basis of vectors that constitute the initial
235:    space, that is, the subspace from which the solver starts to iterate.

237:    Collective on SVD and Vec

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

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

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

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

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

257:    Level: intermediate
258: @*/
259: PetscErrorCode SVDSetInitialSpace(SVD svd,PetscInt n,Vec *is)
260: {
262:   PetscInt       i;
263: 
266:   if (n<0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Argument n cannot be negative");

268:   /* free previous non-processed vectors */
269:   if (svd->nini<0) {
270:     for (i=0;i<-svd->nini;i++) {
271:       VecDestroy(svd->IS[i]);
272:     }
273:     PetscFree(svd->IS);
274:   }

276:   /* get references of passed vectors */
277:   PetscMalloc(n*sizeof(Vec),&svd->IS);
278:   for (i=0;i<n;i++) {
279:     PetscObjectReference((PetscObject)is[i]);
280:     svd->IS[i] = is[i];
281:   }

283:   svd->nini = -n;
284:   svd->setupcalled = 0;
285:   return(0);
286: }