Actual source code: test2.c
slepc-3.14.0 2020-09-30
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2020, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
7: SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9: */
11: static char help[] = "Test dense LME functions.\n\n";
13: #include <slepclme.h>
15: int main(int argc,char **argv)
16: {
18: LME lme;
19: Mat A,B,C,X;
20: PetscInt i,j,n=10,k=2;
21: PetscScalar *As,*Bs,*Cs,*Xs;
22: PetscViewer viewer;
23: PetscBool verbose;
25: SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
26: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
27: PetscOptionsGetInt(NULL,NULL,"-k",&k,NULL);
28: PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
29: PetscPrintf(PETSC_COMM_WORLD,"Dense matrix equations, n=%D.\n",n);
31: /* Create LME object */
32: LMECreate(PETSC_COMM_WORLD,&lme);
34: /* Set up viewer */
35: PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
36: if (verbose) {
37: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);
38: }
40: /* Create matrices */
41: MatCreateSeqDense(PETSC_COMM_SELF,n,n,NULL,&A);
42: PetscObjectSetName((PetscObject)A,"A");
43: MatCreateSeqDense(PETSC_COMM_SELF,n,n,NULL,&B);
44: PetscObjectSetName((PetscObject)B,"B");
45: MatCreateSeqDense(PETSC_COMM_SELF,n,k,NULL,&C);
46: PetscObjectSetName((PetscObject)C,"C");
47: MatCreateSeqDense(PETSC_COMM_SELF,n,n,NULL,&X);
48: PetscObjectSetName((PetscObject)X,"X");
50: /* Fill A with an upper Hessenberg Toeplitz matrix */
51: MatDenseGetArray(A,&As);
52: for (i=0;i<n;i++) As[i+i*n]=3.0-(PetscReal)n/2;
53: for (i=0;i<n-1;i++) As[i+1+i*n]=0.5;
54: for (j=1;j<3;j++) {
55: for (i=0;i<n-j;i++) As[i+(i+j)*n]=1.0;
56: }
57: MatDenseRestoreArray(A,&As);
59: if (verbose) {
60: PetscPrintf(PETSC_COMM_WORLD,"Matrix A - - - - - - - -\n");
61: MatView(A,viewer);
62: }
64: /* Fill B with the 1-D Laplacian matrix */
65: MatDenseGetArray(B,&Bs);
66: for (i=0;i<n;i++) Bs[i+i*n]=2.0;
67: for (i=0;i<n-1;i++) { Bs[i+1+i*n]=-1; Bs[i+(i+1)*n]=-1; }
68: MatDenseRestoreArray(B,&Bs);
70: if (verbose) {
71: PetscPrintf(PETSC_COMM_WORLD,"Matrix B - - - - - - - -\n");
72: MatView(B,viewer);
73: }
75: /* Solve Lyapunov equation A*X+X*A'= -B */
76: PetscPrintf(PETSC_COMM_WORLD,"Solving Lyapunov equation for B\n");
77: MatDenseGetArray(A,&As);
78: MatDenseGetArray(B,&Bs);
79: MatDenseGetArray(X,&Xs);
80: LMEDenseLyapunov(lme,n,As,n,Bs,n,Xs,n);
81: MatDenseRestoreArray(A,&As);
82: MatDenseRestoreArray(B,&Bs);
83: MatDenseRestoreArray(X,&Xs);
84: if (verbose) {
85: PetscPrintf(PETSC_COMM_WORLD,"Solution X - - - - - - - -\n");
86: MatView(X,viewer);
87: }
89: /* Fill C with a full-rank nx2 matrix */
90: MatDenseGetArray(C,&Cs);
91: for (i=0;i<k;i++) Cs[i+i*n] = (i%2)? -1.0: 1.0;
92: MatDenseRestoreArray(C,&Cs);
94: if (verbose) {
95: PetscPrintf(PETSC_COMM_WORLD,"Matrix C - - - - - - - -\n");
96: MatView(C,viewer);
97: }
99: /* Solve Lyapunov equation A*X+X*A'= -C*C' */
100: PetscPrintf(PETSC_COMM_WORLD,"Solving Lyapunov equation for C (Cholesky)\n");
101: MatDenseGetArray(A,&As);
102: MatDenseGetArray(C,&Cs);
103: MatDenseGetArray(X,&Xs);
104: LMEDenseHessLyapunovChol(lme,n,As,n,2,Cs,n,Xs,n,NULL);
105: MatDenseRestoreArray(A,&As);
106: MatDenseRestoreArray(C,&Cs);
107: MatDenseRestoreArray(X,&Xs);
108: if (verbose) {
109: PetscPrintf(PETSC_COMM_WORLD,"Solution X - - - - - - - -\n");
110: MatView(X,viewer);
111: }
113: MatDestroy(&A);
114: MatDestroy(&B);
115: MatDestroy(&C);
116: MatDestroy(&X);
117: LMEDestroy(&lme);
118: SlepcFinalize();
119: return ierr;
120: }
122: /*TEST
124: test:
125: args: -info :lme
126: requires: double
127: filter: sed -e "s/[0-9]\.[0-9]*e[+-]\([0-9]*\)/1e-\\1/g"
129: TEST*/