Actual source code: eige.c
2: #include src/ksp/ksp/kspimpl.h
6: /*@
7: KSPComputeExplicitOperator - Computes the explicit preconditioned operator.
9: Collective on KSP
11: Input Parameter:
12: . ksp - the Krylov subspace context
14: Output Parameter:
15: . mat - the explict preconditioned operator
17: Notes:
18: This computation is done by applying the operators to columns of the
19: identity matrix.
21: Currently, this routine uses a dense matrix format when 1 processor
22: is used and a sparse format otherwise. This routine is costly in general,
23: and is recommended for use only with relatively small systems.
25: Level: advanced
26:
27: .keywords: KSP, compute, explicit, operator
29: .seealso: KSPComputeEigenvaluesExplicitly()
30: @*/
31: PetscErrorCode KSPComputeExplicitOperator(KSP ksp,Mat *mat)
32: {
33: Vec in,out;
35: PetscMPIInt size;
36: PetscInt i,M,m,*rows,start,end;
37: Mat A;
38: MPI_Comm comm;
39: PetscScalar *array,zero = 0.0,one = 1.0;
44: comm = ksp->comm;
46: MPI_Comm_size(comm,&size);
48: VecDuplicate(ksp->vec_sol,&in);
49: VecDuplicate(ksp->vec_sol,&out);
50: VecGetSize(in,&M);
51: VecGetLocalSize(in,&m);
52: VecGetOwnershipRange(in,&start,&end);
53: PetscMalloc((m+1)*sizeof(PetscInt),&rows);
54: for (i=0; i<m; i++) {rows[i] = start + i;}
56: MatCreate(comm,m,m,M,M,mat);
57: if (size == 1) {
58: MatSetType(*mat,MATSEQDENSE);
59: MatSeqDenseSetPreallocation(*mat,PETSC_NULL);
60: } else {
61: MatSetType(*mat,MATMPIAIJ);
62: MatMPIAIJSetPreallocation(*mat,0,PETSC_NULL,0,PETSC_NULL);
63: }
64:
65: PCGetOperators(ksp->pc,&A,PETSC_NULL,PETSC_NULL);
67: for (i=0; i<M; i++) {
69: VecSet(&zero,in);
70: VecSetValues(in,1,&i,&one,INSERT_VALUES);
71: VecAssemblyBegin(in);
72: VecAssemblyEnd(in);
74: KSP_MatMult(ksp,A,in,out);
75: KSP_PCApply(ksp,out,in);
76:
77: VecGetArray(in,&array);
78: MatSetValues(*mat,m,rows,1,&i,array,INSERT_VALUES);
79: VecRestoreArray(in,&array);
81: }
82: PetscFree(rows);
83: VecDestroy(in);
84: VecDestroy(out);
85: MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);
86: MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);
87: return(0);
88: }
90: #include petscblaslapack.h
94: /*@
95: KSPComputeEigenvaluesExplicitly - Computes all of the eigenvalues of the
96: preconditioned operator using LAPACK.
98: Collective on KSP
100: Input Parameter:
101: + ksp - iterative context obtained from KSPCreate()
102: - n - size of arrays r and c
104: Output Parameters:
105: + r - real part of computed eigenvalues
106: - c - complex part of computed eigenvalues
108: Notes:
109: This approach is very slow but will generally provide accurate eigenvalue
110: estimates. This routine explicitly forms a dense matrix representing
111: the preconditioned operator, and thus will run only for relatively small
112: problems, say n < 500.
114: Many users may just want to use the monitoring routine
115: KSPSingularValueMonitor() (which can be set with option -ksp_singmonitor)
116: to print the singular values at each iteration of the linear solve.
118: The preconditoner operator, rhs vector, solution vectors should be
119: set before this routine is called. i.e use KSPSetOperators(),KSPSolve() or
120: KSPSetOperators()
122: Level: advanced
124: .keywords: KSP, compute, eigenvalues, explicitly
126: .seealso: KSPComputeEigenvalues(), KSPSingularValueMonitor(), KSPComputeExtremeSingularValues(), KSPSetOperators(), KSPSolve()
127: @*/
128: PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP ksp,PetscInt nmax,PetscReal *r,PetscReal *c)
129: {
130: Mat BA;
131: PetscErrorCode ierr;
132: PetscBLASInt lierr;
133: PetscMPIInt size,rank;
134: MPI_Comm comm = ksp->comm;
135: PetscScalar *array;
136: Mat A;
137: PetscInt m,row,nz,i,n,dummy;
138: const PetscInt *cols;
139: const PetscScalar *vals;
142: KSPComputeExplicitOperator(ksp,&BA);
143: MPI_Comm_size(comm,&size);
144: MPI_Comm_rank(comm,&rank);
146: MatGetSize(BA,&n,&n);
147: if (size > 1) { /* assemble matrix on first processor */
148: if (!rank) {
149: MatCreate(ksp->comm,n,n,n,n,&A);
150: } else {
151: MatCreate(ksp->comm,0,n,n,n,&A);
152: }
153: MatSetType(A,MATMPIDENSE);
154: MatMPIDenseSetPreallocation(A,PETSC_NULL);
155: PetscLogObjectParent(BA,A);
157: MatGetOwnershipRange(BA,&row,&dummy);
158: MatGetLocalSize(BA,&m,&dummy);
159: for (i=0; i<m; i++) {
160: MatGetRow(BA,row,&nz,&cols,&vals);
161: MatSetValues(A,1,&row,nz,cols,vals,INSERT_VALUES);
162: MatRestoreRow(BA,row,&nz,&cols,&vals);
163: row++;
164: }
166: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
167: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
168: MatGetArray(A,&array);
169: } else {
170: MatGetArray(BA,&array);
171: }
173: #if defined(PETSC_HAVE_ESSL)
174: /* ESSL has a different calling sequence for dgeev() and zgeev() than standard LAPACK */
175: if (!rank) {
176: PetscScalar sdummy,*cwork;
177: PetscReal *work,*realpart;
178: PetscInt clen,idummy,lwork,*perm,zero;
180: #if !defined(PETSC_USE_COMPLEX)
181: clen = n;
182: #else
183: clen = 2*n;
184: #endif
185: PetscMalloc(clen*sizeof(PetscScalar),&cwork);
186: idummy = n;
187: lwork = 5*n;
188: PetscMalloc(lwork*sizeof(PetscReal),&work);
189: PetscMalloc(n*sizeof(PetscReal),&realpart);
190: zero = 0;
191: LAgeev_(&zero,array,&n,cwork,&sdummy,&idummy,&idummy,&n,work,&lwork);
192: PetscFree(work);
194: /* For now we stick with the convention of storing the real and imaginary
195: components of evalues separately. But is this what we really want? */
196: PetscMalloc(n*sizeof(PetscInt),&perm);
198: #if !defined(PETSC_USE_COMPLEX)
199: for (i=0; i<n; i++) {
200: realpart[i] = cwork[2*i];
201: perm[i] = i;
202: }
203: PetscSortRealWithPermutation(n,realpart,perm);
204: for (i=0; i<n; i++) {
205: r[i] = cwork[2*perm[i]];
206: c[i] = cwork[2*perm[i]+1];
207: }
208: #else
209: for (i=0; i<n; i++) {
210: realpart[i] = PetscRealPart(cwork[i]);
211: perm[i] = i;
212: }
213: PetscSortRealWithPermutation(n,realpart,perm);
214: for (i=0; i<n; i++) {
215: r[i] = PetscRealPart(cwork[perm[i]]);
216: c[i] = PetscImaginaryPart(cwork[perm[i]]);
217: }
218: #endif
219: PetscFree(perm);
220: PetscFree(realpart);
221: PetscFree(cwork);
222: }
223: #elif !defined(PETSC_USE_COMPLEX)
224: if (!rank) {
225: PetscScalar *work;
226: PetscReal *realpart,*imagpart;
227: PetscBLASInt idummy,lwork;
228: PetscInt *perm;
230: idummy = n;
231: lwork = 5*n;
232: PetscMalloc(2*n*sizeof(PetscReal),&realpart);
233: imagpart = realpart + n;
234: PetscMalloc(5*n*sizeof(PetscReal),&work);
235: #if defined(PETSC_MISSING_LAPACK_GEEV)
236: SETERRQ(PETSC_ERR_SUP,"GEEV - Lapack routine is unavailable\nNot able to provide eigen values.");
237: #else
238: {
239: PetscBLASInt bn = (PetscBLASInt) n;
240: PetscScalar sdummy;
241: LAgeev_("N","N",&bn,array,&bn,realpart,imagpart,&sdummy,&idummy,&sdummy,&idummy,work,&lwork,&lierr);
242: }
243: #endif
244: if (lierr) SETERRQ1(PETSC_ERR_LIB,"Error in LAPACK routine %d",(int)lierr);
245: PetscFree(work);
246: PetscMalloc(n*sizeof(PetscInt),&perm);
247: for (i=0; i<n; i++) { perm[i] = i;}
248: PetscSortRealWithPermutation(n,realpart,perm);
249: for (i=0; i<n; i++) {
250: r[i] = realpart[perm[i]];
251: c[i] = imagpart[perm[i]];
252: }
253: PetscFree(perm);
254: PetscFree(realpart);
255: }
256: #else
257: if (!rank) {
258: PetscScalar *work,*eigs;
259: PetscReal *rwork;
260: PetscBLASInt idummy,lwork;
261: PetscInt *perm;
263: idummy = n;
264: lwork = 5*n;
265: PetscMalloc(5*n*sizeof(PetscScalar),&work);
266: PetscMalloc(2*n*sizeof(PetscReal),&rwork);
267: PetscMalloc(n*sizeof(PetscScalar),&eigs);
268: #if defined(PETSC_MISSING_LAPACK_GEEV)
269: SETERRQ(PETSC_ERR_SUP,"GEEV - Lapack routine is unavailable\nNot able to provide eigen values.");
270: #else
271: {
272: PetscScalar sdummy;
273: LAgeev_("N","N",&n,array,&n,eigs,&sdummy,&idummy,&sdummy,&idummy,work,&lwork,rwork,&lierr);
274: }
275: #endif
276: if (lierr) SETERRQ1(PETSC_ERR_LIB,"Error in LAPACK routine %d",(int)lierr);
277: PetscFree(work);
278: PetscFree(rwork);
279: PetscMalloc(n*sizeof(PetscInt),&perm);
280: for (i=0; i<n; i++) { perm[i] = i;}
281: for (i=0; i<n; i++) { r[i] = PetscRealPart(eigs[i]);}
282: PetscSortRealWithPermutation(n,r,perm);
283: for (i=0; i<n; i++) {
284: r[i] = PetscRealPart(eigs[perm[i]]);
285: c[i] = PetscImaginaryPart(eigs[perm[i]]);
286: }
287: PetscFree(perm);
288: PetscFree(eigs);
289: }
290: #endif
291: if (size > 1) {
292: MatRestoreArray(A,&array);
293: MatDestroy(A);
294: } else {
295: MatRestoreArray(BA,&array);
296: }
297: MatDestroy(BA);
298: return(0);
299: }