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