Actual source code: cross.c
1: /*
3: SLEPc singular value solver: "cross"
5: Method: Uses a Hermitian eigensolver for A^T*A
7: Last update: Jun 2007
9: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
10: SLEPc - Scalable Library for Eigenvalue Problem Computations
11: Copyright (c) 2002-2009, Universidad Politecnica de Valencia, Spain
13: This file is part of SLEPc.
14:
15: SLEPc is free software: you can redistribute it and/or modify it under the
16: terms of version 3 of the GNU Lesser General Public License as published by
17: the Free Software Foundation.
19: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
20: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
21: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
22: more details.
24: You should have received a copy of the GNU Lesser General Public License
25: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
26: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
27: */
29: #include private/svdimpl.h
30: #include slepceps.h
32: typedef struct {
33: EPS eps;
34: Mat mat;
35: Vec w,diag;
36: } SVD_CROSS;
40: PetscErrorCode ShellMatMult_CROSS(Mat B,Vec x, Vec y)
41: {
43: SVD svd;
44: SVD_CROSS *cross;
45:
47: MatShellGetContext(B,(void**)&svd);
48: cross = (SVD_CROSS *)svd->data;
49: SVDMatMult(svd,PETSC_FALSE,x,cross->w);
50: SVDMatMult(svd,PETSC_TRUE,cross->w,y);
51: return(0);
52: }
56: PetscErrorCode ShellMatGetDiagonal_CROSS(Mat B,Vec d)
57: {
58: PetscErrorCode ierr;
59: SVD svd;
60: SVD_CROSS *cross;
61: PetscInt N,n,i,j,start,end,ncols;
62: PetscScalar *work1,*work2,*diag;
63: const PetscInt *cols;
64: const PetscScalar *vals;
65:
67: MatShellGetContext(B,(void**)&svd);
68: cross = (SVD_CROSS *)svd->data;
69: if (!cross->diag) {
70: /* compute diagonal from rows and store in cross->diag */
71: VecDuplicate(d,&cross->diag);
72: SVDMatGetSize(svd,PETSC_NULL,&N);
73: SVDMatGetLocalSize(svd,PETSC_NULL,&n);
74: PetscMalloc(sizeof(PetscScalar)*N,&work1);
75: PetscMalloc(sizeof(PetscScalar)*N,&work2);
76: for (i=0;i<n;i++) work1[i] = work2[i] = 0.0;
77: if (svd->AT) {
78: MatGetOwnershipRange(svd->AT,&start,&end);
79: for (i=start;i<end;i++) {
80: MatGetRow(svd->AT,i,&ncols,PETSC_NULL,&vals);
81: for (j=0;j<ncols;j++)
82: work1[i] += vals[j]*vals[j];
83: MatRestoreRow(svd->AT,i,&ncols,PETSC_NULL,&vals);
84: }
85: } else {
86: MatGetOwnershipRange(svd->A,&start,&end);
87: for (i=start;i<end;i++) {
88: MatGetRow(svd->A,i,&ncols,&cols,&vals);
89: for (j=0;j<ncols;j++)
90: work1[cols[j]] += vals[j]*vals[j];
91: MatRestoreRow(svd->A,i,&ncols,&cols,&vals);
92: }
93: }
94: MPI_Allreduce(work1,work2,N,MPIU_SCALAR,MPI_SUM,((PetscObject)svd)->comm);
95: VecGetOwnershipRange(cross->diag,&start,&end);
96: VecGetArray(cross->diag,&diag);
97: for (i=start;i<end;i++)
98: diag[i-start] = work2[i];
99: VecRestoreArray(cross->diag,&diag);
100: PetscFree(work1);
101: PetscFree(work2);
102: }
103: VecCopy(cross->diag,d);
104: return(0);
105: }
109: PetscErrorCode SVDSetUp_CROSS(SVD svd)
110: {
111: PetscErrorCode ierr;
112: SVD_CROSS *cross = (SVD_CROSS *)svd->data;
113: PetscInt n;
116: if (cross->mat) {
117: MatDestroy(cross->mat);
118: VecDestroy(cross->w);
119: }
120: if (cross->diag) {
121: VecDestroy(cross->diag);
122: }
123:
124: SVDMatGetLocalSize(svd,PETSC_NULL,&n);
125: MatCreateShell(((PetscObject)svd)->comm,n,n,PETSC_DETERMINE,PETSC_DETERMINE,svd,&cross->mat);
126: MatShellSetOperation(cross->mat,MATOP_MULT,(void(*)(void))ShellMatMult_CROSS);
127: MatShellSetOperation(cross->mat,MATOP_GET_DIAGONAL,(void(*)(void))ShellMatGetDiagonal_CROSS);
128: SVDMatGetVecs(svd,PETSC_NULL,&cross->w);
130: EPSSetOperators(cross->eps,cross->mat,PETSC_NULL);
131: EPSSetProblemType(cross->eps,EPS_HEP);
132: EPSSetWhichEigenpairs(cross->eps,svd->which == SVD_LARGEST ? EPS_LARGEST_REAL : EPS_SMALLEST_REAL);
133: EPSSetDimensions(cross->eps,svd->nsv,svd->ncv,svd->mpd);
134: EPSSetTolerances(cross->eps,svd->tol,svd->max_it);
135: EPSSetUp(cross->eps);
136: EPSGetDimensions(cross->eps,PETSC_NULL,&svd->ncv,&svd->mpd);
137: EPSGetTolerances(cross->eps,&svd->tol,&svd->max_it);
138: return(0);
139: }
143: PetscErrorCode SVDSolve_CROSS(SVD svd)
144: {
146: SVD_CROSS *cross = (SVD_CROSS *)svd->data;
147: PetscInt i;
148: PetscScalar sigma;
149:
151: EPSSetInitialVector(cross->eps,svd->vec_initial);
152: EPSSolve(cross->eps);
153: EPSGetConverged(cross->eps,&svd->nconv);
154: EPSGetIterationNumber(cross->eps,&svd->its);
155: EPSGetConvergedReason(cross->eps,(EPSConvergedReason*)&svd->reason);
156: for (i=0;i<svd->nconv;i++) {
157: EPSGetEigenpair(cross->eps,i,&sigma,PETSC_NULL,svd->V[i],PETSC_NULL);
158: svd->sigma[i] = sqrt(PetscRealPart(sigma));
159: }
160: return(0);
161: }
165: PetscErrorCode SVDMonitor_CROSS(EPS eps,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *ctx)
166: {
167: PetscInt i;
168: SVD svd = (SVD)ctx;
171: for (i=0;i<nest;i++) {
172: svd->sigma[i] = sqrt(PetscRealPart(eigr[i]));
173: svd->errest[i] = errest[i];
174: }
175: SVDMonitor(svd,its,nconv,svd->sigma,svd->errest,nest);
176: return(0);
177: }
181: PetscErrorCode SVDSetFromOptions_CROSS(SVD svd)
182: {
184: SVD_CROSS *cross = (SVD_CROSS *)svd->data;
187: EPSSetFromOptions(cross->eps);
188: return(0);
189: }
194: PetscErrorCode SVDCrossSetEPS_CROSS(SVD svd,EPS eps)
195: {
197: SVD_CROSS *cross = (SVD_CROSS *)svd->data;
202: PetscObjectReference((PetscObject)eps);
203: EPSDestroy(cross->eps);
204: cross->eps = eps;
205: svd->setupcalled = 0;
206: return(0);
207: }
212: /*@
213: SVDCrossSetEPS - Associate an eigensolver object (EPS) to the
214: singular value solver.
216: Collective on SVD
218: Input Parameters:
219: + svd - singular value solver
220: - eps - the eigensolver object
222: Level: advanced
224: .seealso: SVDCrossGetEPS()
225: @*/
226: PetscErrorCode SVDCrossSetEPS(SVD svd,EPS eps)
227: {
228: PetscErrorCode ierr, (*f)(SVD,EPS eps);
232: PetscObjectQueryFunction((PetscObject)svd,"SVDCrossSetEPS_C",(void (**)())&f);
233: if (f) {
234: (*f)(svd,eps);
235: }
236: return(0);
237: }
242: PetscErrorCode SVDCrossGetEPS_CROSS(SVD svd,EPS *eps)
243: {
244: SVD_CROSS *cross = (SVD_CROSS *)svd->data;
248: *eps = cross->eps;
249: return(0);
250: }
255: /*@C
256: SVDCrossGetEPS - Retrieve the eigensolver object (EPS) associated
257: to the singular value solver.
259: Not Collective
261: Input Parameter:
262: . svd - singular value solver
264: Output Parameter:
265: . eps - the eigensolver object
267: Level: advanced
269: .seealso: SVDCrossSetEPS()
270: @*/
271: PetscErrorCode SVDCrossGetEPS(SVD svd,EPS *eps)
272: {
273: PetscErrorCode ierr, (*f)(SVD,EPS *eps);
277: PetscObjectQueryFunction((PetscObject)svd,"SVDCrossGetEPS_C",(void (**)())&f);
278: if (f) {
279: (*f)(svd,eps);
280: }
281: return(0);
282: }
286: PetscErrorCode SVDView_CROSS(SVD svd,PetscViewer viewer)
287: {
289: SVD_CROSS *cross = (SVD_CROSS *)svd->data;
292: EPSView(cross->eps,viewer);
293: return(0);
294: }
298: PetscErrorCode SVDDestroy_CROSS(SVD svd)
299: {
301: SVD_CROSS *cross = (SVD_CROSS *)svd->data;
304: EPSDestroy(cross->eps);
305: if (cross->mat) {
306: MatDestroy(cross->mat);
307: VecDestroy(cross->w);
308: }
309: if (cross->diag) {
310: VecDestroy(cross->diag);
311: }
312: PetscFree(svd->data);
313: return(0);
314: }
319: PetscErrorCode SVDCreate_CROSS(SVD svd)
320: {
322: SVD_CROSS *cross;
323: ST st;
324:
326: PetscNew(SVD_CROSS,&cross);
327: PetscLogObjectMemory(svd,sizeof(SVD_CROSS));
328: svd->data = (void *)cross;
329: svd->ops->solve = SVDSolve_CROSS;
330: svd->ops->setup = SVDSetUp_CROSS;
331: svd->ops->setfromoptions = SVDSetFromOptions_CROSS;
332: svd->ops->destroy = SVDDestroy_CROSS;
333: svd->ops->view = SVDView_CROSS;
334: PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDCrossSetEPS_C","SVDCrossSetEPS_CROSS",SVDCrossSetEPS_CROSS);
335: PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDCrossGetEPS_C","SVDCrossGetEPS_CROSS",SVDCrossGetEPS_CROSS);
337: EPSCreate(((PetscObject)svd)->comm,&cross->eps);
338: EPSSetOptionsPrefix(cross->eps,((PetscObject)svd)->prefix);
339: EPSAppendOptionsPrefix(cross->eps,"svd_");
340: PetscObjectIncrementTabLevel((PetscObject)cross->eps,(PetscObject)svd,1);
341: PetscLogObjectParent(svd,cross->eps);
342: EPSSetIP(cross->eps,svd->ip);
343: EPSSetWhichEigenpairs(cross->eps,EPS_LARGEST_REAL);
344: EPSMonitorSet(cross->eps,SVDMonitor_CROSS,svd,PETSC_NULL);
345: EPSGetST(cross->eps,&st);
346: STSetMatMode(st,STMATMODE_SHELL);
347: cross->mat = PETSC_NULL;
348: cross->w = PETSC_NULL;
349: cross->diag = PETSC_NULL;
350: return(0);
351: }