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-2007, Universidad Politecnica de Valencia, Spain
13: This file is part of SLEPc. See the README file for conditions of use
14: and additional information.
15: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
16: */
18: #include src/svd/svdimpl.h
19: #include slepceps.h
21: typedef struct {
22: EPS eps;
23: Mat mat;
24: Vec w,diag;
25: } SVD_CROSS;
29: PetscErrorCode ShellMatMult_CROSS(Mat B,Vec x, Vec y)
30: {
32: SVD svd;
33: SVD_CROSS *cross;
34:
36: MatShellGetContext(B,(void**)&svd);
37: cross = (SVD_CROSS *)svd->data;
38: SVDMatMult(svd,PETSC_FALSE,x,cross->w);
39: SVDMatMult(svd,PETSC_TRUE,cross->w,y);
40: return(0);
41: }
45: PetscErrorCode ShellMatGetDiagonal_CROSS(Mat B,Vec d)
46: {
47: PetscErrorCode ierr;
48: SVD svd;
49: SVD_CROSS *cross;
50: PetscInt N,n,i,j,start,end,ncols;
51: PetscScalar *work1,*work2,*diag;
52: const PetscInt *cols;
53: const PetscScalar *vals;
54:
56: MatShellGetContext(B,(void**)&svd);
57: cross = (SVD_CROSS *)svd->data;
58: if (!cross->diag) {
59: /* compute diagonal from rows and store in cross->diag */
60: VecDuplicate(d,&cross->diag);
61: SVDMatGetSize(svd,PETSC_NULL,&N);
62: SVDMatGetLocalSize(svd,PETSC_NULL,&n);
63: PetscMalloc(sizeof(PetscScalar)*N,&work1);
64: PetscMalloc(sizeof(PetscScalar)*N,&work2);
65: for (i=0;i<n;i++) work1[i] = work2[i] = 0.0;
66: if (svd->AT) {
67: MatGetOwnershipRange(svd->AT,&start,&end);
68: for (i=start;i<end;i++) {
69: MatGetRow(svd->AT,i,&ncols,PETSC_NULL,&vals);
70: for (j=0;j<ncols;j++)
71: work1[i] += vals[j]*vals[j];
72: MatRestoreRow(svd->AT,i,&ncols,PETSC_NULL,&vals);
73: }
74: } else {
75: MatGetOwnershipRange(svd->A,&start,&end);
76: for (i=start;i<end;i++) {
77: MatGetRow(svd->A,i,&ncols,&cols,&vals);
78: for (j=0;j<ncols;j++)
79: work1[cols[j]] += vals[j]*vals[j];
80: MatRestoreRow(svd->A,i,&ncols,&cols,&vals);
81: }
82: }
83: MPI_Allreduce(work1,work2,N,MPIU_SCALAR,MPI_SUM,svd->comm);
84: VecGetOwnershipRange(cross->diag,&start,&end);
85: VecGetArray(cross->diag,&diag);
86: for (i=start;i<end;i++)
87: diag[i-start] = work2[i];
88: VecRestoreArray(cross->diag,&diag);
89: PetscFree(work1);
90: PetscFree(work2);
91: }
92: VecCopy(cross->diag,d);
93: return(0);
94: }
98: PetscErrorCode SVDSetUp_CROSS(SVD svd)
99: {
100: PetscErrorCode ierr;
101: SVD_CROSS *cross = (SVD_CROSS *)svd->data;
102: PetscInt n;
105: if (cross->mat) {
106: MatDestroy(cross->mat);
107: VecDestroy(cross->w);
108: }
109: if (cross->diag) {
110: VecDestroy(cross->diag);
111: }
112:
113: SVDMatGetLocalSize(svd,PETSC_NULL,&n);
114: MatCreateShell(svd->comm,n,n,PETSC_DETERMINE,PETSC_DETERMINE,svd,&cross->mat);
115: MatShellSetOperation(cross->mat,MATOP_MULT,(void(*)(void))ShellMatMult_CROSS);
116: MatShellSetOperation(cross->mat,MATOP_GET_DIAGONAL,(void(*)(void))ShellMatGetDiagonal_CROSS);
117: SVDMatGetVecs(svd,PETSC_NULL,&cross->w);
119: EPSSetOperators(cross->eps,cross->mat,PETSC_NULL);
120: EPSSetProblemType(cross->eps,EPS_HEP);
121: EPSSetDimensions(cross->eps,svd->nsv,svd->ncv);
122: EPSSetTolerances(cross->eps,svd->tol,svd->max_it);
123: EPSSetUp(cross->eps);
124: EPSGetDimensions(cross->eps,PETSC_NULL,&svd->ncv);
125: EPSGetTolerances(cross->eps,&svd->tol,&svd->max_it);
126: return(0);
127: }
131: PetscErrorCode SVDSolve_CROSS(SVD svd)
132: {
134: SVD_CROSS *cross = (SVD_CROSS *)svd->data;
135: int i;
136: PetscScalar sigma;
137:
139: EPSSetWhichEigenpairs(cross->eps,svd->which == SVD_LARGEST ? EPS_LARGEST_REAL : EPS_SMALLEST_REAL);
140: EPSSetInitialVector(cross->eps,svd->vec_initial);
141: EPSSolve(cross->eps);
142: EPSGetConverged(cross->eps,&svd->nconv);
143: EPSGetIterationNumber(cross->eps,&svd->its);
144: EPSGetConvergedReason(cross->eps,(EPSConvergedReason*)&svd->reason);
145: for (i=0;i<svd->nconv;i++) {
146: EPSGetEigenpair(cross->eps,i,&sigma,PETSC_NULL,svd->V[i],PETSC_NULL);
147: svd->sigma[i] = sqrt(PetscRealPart(sigma));
148: }
149: return(0);
150: }
154: PetscErrorCode SVDMonitor_CROSS(EPS eps,int its,int nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,int nest,void *ctx)
155: {
156: int i;
157: SVD svd = (SVD)ctx;
160: for (i=0;i<nest;i++) {
161: svd->sigma[i] = sqrt(PetscRealPart(eigr[i]));
162: svd->errest[i] = errest[i];
163: }
164: SVDMonitor(svd,its,nconv,svd->sigma,svd->errest,nest);
165: return(0);
166: }
170: PetscErrorCode SVDSetFromOptions_CROSS(SVD svd)
171: {
173: SVD_CROSS *cross = (SVD_CROSS *)svd->data;
176: EPSSetFromOptions(cross->eps);
177: return(0);
178: }
183: PetscErrorCode SVDCrossSetEPS_CROSS(SVD svd,EPS eps)
184: {
186: SVD_CROSS *cross = (SVD_CROSS *)svd->data;
191: PetscObjectReference((PetscObject)eps);
192: EPSDestroy(cross->eps);
193: cross->eps = eps;
194: svd->setupcalled = 0;
195: return(0);
196: }
201: /*@
202: SVDCrossSetEPS - Associate an eigensolver object (EPS) to the
203: singular value solver.
205: Collective on SVD
207: Input Parameters:
208: + svd - singular value solver
209: - eps - the eigensolver object
211: Level: advanced
213: .seealso: SVDCrossGetEPS()
214: @*/
215: PetscErrorCode SVDCrossSetEPS(SVD svd,EPS eps)
216: {
217: PetscErrorCode ierr, (*f)(SVD,EPS eps);
221: PetscObjectQueryFunction((PetscObject)svd,"SVDCrossSetEPS_C",(void (**)())&f);
222: if (f) {
223: (*f)(svd,eps);
224: }
225: return(0);
226: }
231: PetscErrorCode SVDCrossGetEPS_CROSS(SVD svd,EPS *eps)
232: {
233: SVD_CROSS *cross = (SVD_CROSS *)svd->data;
237: *eps = cross->eps;
238: return(0);
239: }
244: /*@C
245: SVDCrossGetEPS - Retrieve the eigensolver object (EPS) associated
246: to the singular value solver.
248: Not Collective
250: Input Parameter:
251: . svd - singular value solver
253: Output Parameter:
254: . eps - the eigensolver object
256: Level: advanced
258: .seealso: SVDCrossSetEPS()
259: @*/
260: PetscErrorCode SVDCrossGetEPS(SVD svd,EPS *eps)
261: {
262: PetscErrorCode ierr, (*f)(SVD,EPS *eps);
266: PetscObjectQueryFunction((PetscObject)svd,"SVDCrossGetEPS_C",(void (**)())&f);
267: if (f) {
268: (*f)(svd,eps);
269: }
270: return(0);
271: }
275: PetscErrorCode SVDView_CROSS(SVD svd,PetscViewer viewer)
276: {
278: SVD_CROSS *cross = (SVD_CROSS *)svd->data;
281: EPSView(cross->eps,viewer);
282: return(0);
283: }
287: PetscErrorCode SVDDestroy_CROSS(SVD svd)
288: {
290: SVD_CROSS *cross = (SVD_CROSS *)svd->data;
293: EPSDestroy(cross->eps);
294: if (cross->mat) {
295: MatDestroy(cross->mat);
296: VecDestroy(cross->w);
297: }
298: if (cross->diag) {
299: VecDestroy(cross->diag);
300: }
301: PetscFree(svd->data);
302: return(0);
303: }
308: PetscErrorCode SVDCreate_CROSS(SVD svd)
309: {
311: SVD_CROSS *cross;
312: ST st;
313:
315: PetscNew(SVD_CROSS,&cross);
316: PetscLogObjectMemory(svd,sizeof(SVD_CROSS));
317: svd->data = (void *)cross;
318: svd->ops->solve = SVDSolve_CROSS;
319: svd->ops->setup = SVDSetUp_CROSS;
320: svd->ops->setfromoptions = SVDSetFromOptions_CROSS;
321: svd->ops->destroy = SVDDestroy_CROSS;
322: svd->ops->view = SVDView_CROSS;
323: PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDCrossSetEPS_C","SVDCrossSetEPS_CROSS",SVDCrossSetEPS_CROSS);
324: PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDCrossGetEPS_C","SVDCrossGetEPS_CROSS",SVDCrossGetEPS_CROSS);
326: EPSCreate(svd->comm,&cross->eps);
327: EPSSetOptionsPrefix(cross->eps,svd->prefix);
328: EPSAppendOptionsPrefix(cross->eps,"svd_");
329: PetscLogObjectParent(svd,cross->eps);
330: EPSSetIP(cross->eps,svd->ip);
331: EPSSetWhichEigenpairs(cross->eps,EPS_LARGEST_REAL);
332: EPSMonitorSet(cross->eps,SVDMonitor_CROSS,svd,PETSC_NULL);
333: EPSGetST(cross->eps,&st);
334: STSetMatMode(st,STMATMODE_SHELL);
335: cross->mat = PETSC_NULL;
336: cross->w = PETSC_NULL;
337: cross->diag = PETSC_NULL;
338: return(0);
339: }