Actual source code: eige.c
2: #include <petsc/private/kspimpl.h>
3: #include <petscdm.h>
4: #include <petscblaslapack.h>
6: typedef struct {
7: KSP ksp;
8: Vec work;
9: } Mat_KSP;
11: static PetscErrorCode MatCreateVecs_KSP(Mat A, Vec *X, Vec *Y)
12: {
13: Mat_KSP *ctx;
14: Mat M;
16: MatShellGetContext(A, &ctx);
17: KSPGetOperators(ctx->ksp, &M, NULL);
18: MatCreateVecs(M, X, Y);
19: return 0;
20: }
22: static PetscErrorCode MatMult_KSP(Mat A, Vec X, Vec Y)
23: {
24: Mat_KSP *ctx;
26: MatShellGetContext(A, &ctx);
27: KSP_PCApplyBAorAB(ctx->ksp, X, Y, ctx->work);
28: return 0;
29: }
31: /*@
32: KSPComputeOperator - Computes the explicit preconditioned operator, including diagonal scaling and null
33: space removal if applicable.
35: Collective on ksp
37: Input Parameters:
38: + ksp - the Krylov subspace context
39: - mattype - the matrix type to be used
41: Output Parameter:
42: . mat - the explicit preconditioned operator
44: Notes:
45: This computation is done by applying the operators to columns of the
46: identity matrix.
48: Currently, this routine uses a dense matrix format for the output operator if mattype == NULL.
49: This routine is costly in general, and is recommended for use only with relatively small systems.
51: Level: advanced
53: .seealso: `KSPComputeEigenvaluesExplicitly()`, `PCComputeOperator()`, `KSPSetDiagonalScale()`, `KSPSetNullSpace()`, `MatType`
54: @*/
55: PetscErrorCode KSPComputeOperator(KSP ksp, MatType mattype, Mat *mat)
56: {
57: PetscInt N, M, m, n;
58: Mat_KSP ctx;
59: Mat A, Aksp;
63: KSPGetOperators(ksp, &A, NULL);
64: MatGetLocalSize(A, &m, &n);
65: MatGetSize(A, &M, &N);
66: MatCreateShell(PetscObjectComm((PetscObject)ksp), m, n, M, N, &ctx, &Aksp);
67: MatShellSetOperation(Aksp, MATOP_MULT, (void (*)(void))MatMult_KSP);
68: MatShellSetOperation(Aksp, MATOP_CREATE_VECS, (void (*)(void))MatCreateVecs_KSP);
69: ctx.ksp = ksp;
70: MatCreateVecs(A, &ctx.work, NULL);
71: MatComputeOperator(Aksp, mattype, mat);
72: VecDestroy(&ctx.work);
73: MatDestroy(&Aksp);
74: return 0;
75: }
77: /*@
78: KSPComputeEigenvaluesExplicitly - Computes all of the eigenvalues of the
79: preconditioned operator using LAPACK.
81: Collective on ksp
83: Input Parameters:
84: + ksp - iterative context obtained from KSPCreate()
85: - n - size of arrays r and c
87: Output Parameters:
88: + r - real part of computed eigenvalues, provided by user with a dimension at least of n
89: - c - complex part of computed eigenvalues, provided by user with a dimension at least of n
91: Notes:
92: This approach is very slow but will generally provide accurate eigenvalue
93: estimates. This routine explicitly forms a dense matrix representing
94: the preconditioned operator, and thus will run only for relatively small
95: problems, say n < 500.
97: Many users may just want to use the monitoring routine
98: KSPMonitorSingularValue() (which can be set with option -ksp_monitor_singular_value)
99: to print the singular values at each iteration of the linear solve.
101: The preconditoner operator, rhs vector, solution vectors should be
102: set before this routine is called. i.e use KSPSetOperators(),KSPSolve() or
103: KSPSetOperators()
105: Level: advanced
107: .seealso: `KSPComputeEigenvalues()`, `KSPMonitorSingularValue()`, `KSPComputeExtremeSingularValues()`, `KSPSetOperators()`, `KSPSolve()`
108: @*/
109: PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP ksp, PetscInt nmax, PetscReal r[], PetscReal c[])
110: {
111: Mat BA;
112: PetscMPIInt size, rank;
113: MPI_Comm comm;
114: PetscScalar *array;
115: Mat A;
116: PetscInt m, row, nz, i, n, dummy;
117: const PetscInt *cols;
118: const PetscScalar *vals;
120: PetscObjectGetComm((PetscObject)ksp, &comm);
121: KSPComputeOperator(ksp, MATDENSE, &BA);
122: MPI_Comm_size(comm, &size);
123: MPI_Comm_rank(comm, &rank);
125: MatGetSize(BA, &n, &n);
126: if (size > 1) { /* assemble matrix on first processor */
127: MatCreate(PetscObjectComm((PetscObject)ksp), &A);
128: if (rank == 0) {
129: MatSetSizes(A, n, n, n, n);
130: } else {
131: MatSetSizes(A, 0, 0, n, n);
132: }
133: MatSetType(A, MATMPIDENSE);
134: MatMPIDenseSetPreallocation(A, NULL);
136: MatGetOwnershipRange(BA, &row, &dummy);
137: MatGetLocalSize(BA, &m, &dummy);
138: for (i = 0; i < m; i++) {
139: MatGetRow(BA, row, &nz, &cols, &vals);
140: MatSetValues(A, 1, &row, nz, cols, vals, INSERT_VALUES);
141: MatRestoreRow(BA, row, &nz, &cols, &vals);
142: row++;
143: }
145: MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);
146: MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);
147: MatDenseGetArray(A, &array);
148: } else {
149: MatDenseGetArray(BA, &array);
150: }
152: #if !defined(PETSC_USE_COMPLEX)
153: if (rank == 0) {
154: PetscScalar *work;
155: PetscReal *realpart, *imagpart;
156: PetscBLASInt idummy, lwork;
157: PetscInt *perm;
159: idummy = n;
160: lwork = 5 * n;
161: PetscMalloc2(n, &realpart, n, &imagpart);
162: PetscMalloc1(5 * n, &work);
163: {
164: PetscBLASInt lierr;
165: PetscScalar sdummy;
166: PetscBLASInt bn;
168: PetscBLASIntCast(n, &bn);
169: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
170: PetscCallBLAS("LAPACKgeev", LAPACKgeev_("N", "N", &bn, array, &bn, realpart, imagpart, &sdummy, &idummy, &sdummy, &idummy, work, &lwork, &lierr));
172: PetscFPTrapPop();
173: }
174: PetscFree(work);
175: PetscMalloc1(n, &perm);
177: for (i = 0; i < n; i++) perm[i] = i;
178: PetscSortRealWithPermutation(n, realpart, perm);
179: for (i = 0; i < n; i++) {
180: r[i] = realpart[perm[i]];
181: c[i] = imagpart[perm[i]];
182: }
183: PetscFree(perm);
184: PetscFree2(realpart, imagpart);
185: }
186: #else
187: if (rank == 0) {
188: PetscScalar *work, *eigs;
189: PetscReal *rwork;
190: PetscBLASInt idummy, lwork;
191: PetscInt *perm;
193: idummy = n;
194: lwork = 5 * n;
195: PetscMalloc1(5 * n, &work);
196: PetscMalloc1(2 * n, &rwork);
197: PetscMalloc1(n, &eigs);
198: {
199: PetscBLASInt lierr;
200: PetscScalar sdummy;
201: PetscBLASInt nb;
202: PetscBLASIntCast(n, &nb);
203: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
204: PetscCallBLAS("LAPACKgeev", LAPACKgeev_("N", "N", &nb, array, &nb, eigs, &sdummy, &idummy, &sdummy, &idummy, work, &lwork, rwork, &lierr));
206: PetscFPTrapPop();
207: }
208: PetscFree(work);
209: PetscFree(rwork);
210: PetscMalloc1(n, &perm);
211: for (i = 0; i < n; i++) perm[i] = i;
212: for (i = 0; i < n; i++) r[i] = PetscRealPart(eigs[i]);
213: PetscSortRealWithPermutation(n, r, perm);
214: for (i = 0; i < n; i++) {
215: r[i] = PetscRealPart(eigs[perm[i]]);
216: c[i] = PetscImaginaryPart(eigs[perm[i]]);
217: }
218: PetscFree(perm);
219: PetscFree(eigs);
220: }
221: #endif
222: if (size > 1) {
223: MatDenseRestoreArray(A, &array);
224: MatDestroy(&A);
225: } else {
226: MatDenseRestoreArray(BA, &array);
227: }
228: MatDestroy(&BA);
229: return 0;
230: }
232: static PetscErrorCode PolyEval(PetscInt nroots, const PetscReal *r, const PetscReal *c, PetscReal x, PetscReal y, PetscReal *px, PetscReal *py)
233: {
234: PetscInt i;
235: PetscReal rprod = 1, iprod = 0;
237: for (i = 0; i < nroots; i++) {
238: PetscReal rnew = rprod * (x - r[i]) - iprod * (y - c[i]);
239: PetscReal inew = rprod * (y - c[i]) + iprod * (x - r[i]);
240: rprod = rnew;
241: iprod = inew;
242: }
243: *px = rprod;
244: *py = iprod;
245: return 0;
246: }
248: #include <petscdraw.h>
249: /* collective on ksp */
250: PetscErrorCode KSPPlotEigenContours_Private(KSP ksp, PetscInt neig, const PetscReal *r, const PetscReal *c)
251: {
252: PetscReal xmin, xmax, ymin, ymax, *xloc, *yloc, *value, px0, py0, rscale, iscale;
253: PetscInt M, N, i, j;
254: PetscMPIInt rank;
255: PetscViewer viewer;
256: PetscDraw draw;
257: PetscDrawAxis drawaxis;
259: MPI_Comm_rank(PetscObjectComm((PetscObject)ksp), &rank);
260: if (rank) return 0;
261: M = 80;
262: N = 80;
263: xmin = r[0];
264: xmax = r[0];
265: ymin = c[0];
266: ymax = c[0];
267: for (i = 1; i < neig; i++) {
268: xmin = PetscMin(xmin, r[i]);
269: xmax = PetscMax(xmax, r[i]);
270: ymin = PetscMin(ymin, c[i]);
271: ymax = PetscMax(ymax, c[i]);
272: }
273: PetscMalloc3(M, &xloc, N, &yloc, M * N, &value);
274: for (i = 0; i < M; i++) xloc[i] = xmin - 0.1 * (xmax - xmin) + 1.2 * (xmax - xmin) * i / (M - 1);
275: for (i = 0; i < N; i++) yloc[i] = ymin - 0.1 * (ymax - ymin) + 1.2 * (ymax - ymin) * i / (N - 1);
276: PolyEval(neig, r, c, 0, 0, &px0, &py0);
277: rscale = px0 / (PetscSqr(px0) + PetscSqr(py0));
278: iscale = -py0 / (PetscSqr(px0) + PetscSqr(py0));
279: for (j = 0; j < N; j++) {
280: for (i = 0; i < M; i++) {
281: PetscReal px, py, tx, ty, tmod;
282: PolyEval(neig, r, c, xloc[i], yloc[j], &px, &py);
283: tx = px * rscale - py * iscale;
284: ty = py * rscale + px * iscale;
285: tmod = PetscSqr(tx) + PetscSqr(ty); /* modulus of the complex polynomial */
286: if (tmod > 1) tmod = 1.0;
287: if (tmod > 0.5 && tmod < 1) tmod = 0.5;
288: if (tmod > 0.2 && tmod < 0.5) tmod = 0.2;
289: if (tmod > 0.05 && tmod < 0.2) tmod = 0.05;
290: if (tmod < 1e-3) tmod = 1e-3;
291: value[i + j * M] = PetscLogReal(tmod) / PetscLogReal(10.0);
292: }
293: }
294: PetscViewerDrawOpen(PETSC_COMM_SELF, NULL, "Iteratively Computed Eigen-contours", PETSC_DECIDE, PETSC_DECIDE, 450, 450, &viewer);
295: PetscViewerDrawGetDraw(viewer, 0, &draw);
296: PetscDrawTensorContour(draw, M, N, NULL, NULL, value);
297: if (0) {
298: PetscDrawAxisCreate(draw, &drawaxis);
299: PetscDrawAxisSetLimits(drawaxis, xmin, xmax, ymin, ymax);
300: PetscDrawAxisSetLabels(drawaxis, "Eigen-counters", "real", "imag");
301: PetscDrawAxisDraw(drawaxis);
302: PetscDrawAxisDestroy(&drawaxis);
303: }
304: PetscViewerDestroy(&viewer);
305: PetscFree3(xloc, yloc, value);
306: return 0;
307: }