Actual source code: test2.c
slepc-3.17.0 2022-03-31
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-, 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: {
17: LME lme;
18: Mat A,B,C,X;
19: PetscInt i,j,n=10,k=2;
20: PetscScalar *As,*Bs,*Cs,*Xs;
21: PetscViewer viewer;
22: PetscBool verbose;
24: SlepcInitialize(&argc,&argv,(char*)0,help);
25: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
26: PetscOptionsGetInt(NULL,NULL,"-k",&k,NULL);
27: PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
28: PetscPrintf(PETSC_COMM_WORLD,"Dense matrix equations, n=%" PetscInt_FMT ".\n",n);
30: /* Create LME object */
31: LMECreate(PETSC_COMM_WORLD,&lme);
33: /* Set up viewer */
34: PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
35: if (verbose) PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);
37: /* Create matrices */
38: MatCreateSeqDense(PETSC_COMM_SELF,n,n,NULL,&A);
39: PetscObjectSetName((PetscObject)A,"A");
40: MatCreateSeqDense(PETSC_COMM_SELF,n,n,NULL,&B);
41: PetscObjectSetName((PetscObject)B,"B");
42: MatCreateSeqDense(PETSC_COMM_SELF,n,k,NULL,&C);
43: PetscObjectSetName((PetscObject)C,"C");
44: MatCreateSeqDense(PETSC_COMM_SELF,n,n,NULL,&X);
45: PetscObjectSetName((PetscObject)X,"X");
47: /* Fill A with an upper Hessenberg Toeplitz matrix */
48: MatDenseGetArray(A,&As);
49: for (i=0;i<n;i++) As[i+i*n]=3.0-(PetscReal)n/2;
50: for (i=0;i<n-1;i++) As[i+1+i*n]=0.5;
51: for (j=1;j<3;j++) {
52: for (i=0;i<n-j;i++) As[i+(i+j)*n]=1.0;
53: }
54: MatDenseRestoreArray(A,&As);
56: if (verbose) {
57: PetscPrintf(PETSC_COMM_WORLD,"Matrix A - - - - - - - -\n");
58: MatView(A,viewer);
59: }
61: /* Fill B with the 1-D Laplacian matrix */
62: MatDenseGetArray(B,&Bs);
63: for (i=0;i<n;i++) Bs[i+i*n]=2.0;
64: for (i=0;i<n-1;i++) { Bs[i+1+i*n]=-1; Bs[i+(i+1)*n]=-1; }
65: MatDenseRestoreArray(B,&Bs);
67: if (verbose) {
68: PetscPrintf(PETSC_COMM_WORLD,"Matrix B - - - - - - - -\n");
69: MatView(B,viewer);
70: }
72: /* Solve Lyapunov equation A*X+X*A'= -B */
73: PetscPrintf(PETSC_COMM_WORLD,"Solving Lyapunov equation for B\n");
74: MatDenseGetArray(A,&As);
75: MatDenseGetArray(B,&Bs);
76: MatDenseGetArray(X,&Xs);
77: LMEDenseLyapunov(lme,n,As,n,Bs,n,Xs,n);
78: MatDenseRestoreArray(A,&As);
79: MatDenseRestoreArray(B,&Bs);
80: MatDenseRestoreArray(X,&Xs);
81: if (verbose) {
82: PetscPrintf(PETSC_COMM_WORLD,"Solution X - - - - - - - -\n");
83: MatView(X,viewer);
84: }
86: /* Fill C with a full-rank nx2 matrix */
87: MatDenseGetArray(C,&Cs);
88: for (i=0;i<k;i++) Cs[i+i*n] = (i%2)? -1.0: 1.0;
89: MatDenseRestoreArray(C,&Cs);
91: if (verbose) {
92: PetscPrintf(PETSC_COMM_WORLD,"Matrix C - - - - - - - -\n");
93: MatView(C,viewer);
94: }
96: /* Solve Lyapunov equation A*X+X*A'= -C*C' */
97: PetscPrintf(PETSC_COMM_WORLD,"Solving Lyapunov equation for C (Cholesky)\n");
98: MatDenseGetArray(A,&As);
99: MatDenseGetArray(C,&Cs);
100: MatDenseGetArray(X,&Xs);
101: LMEDenseHessLyapunovChol(lme,n,As,n,2,Cs,n,Xs,n,NULL);
102: MatDenseRestoreArray(A,&As);
103: MatDenseRestoreArray(C,&Cs);
104: MatDenseRestoreArray(X,&Xs);
105: if (verbose) {
106: PetscPrintf(PETSC_COMM_WORLD,"Solution X - - - - - - - -\n");
107: MatView(X,viewer);
108: }
110: MatDestroy(&A);
111: MatDestroy(&B);
112: MatDestroy(&C);
113: MatDestroy(&X);
114: LMEDestroy(&lme);
115: SlepcFinalize();
116: return 0;
117: }
119: /*TEST
121: test:
122: args: -info :lme
123: requires: double
124: filter: sed -e "s/[0-9]\.[0-9]*e[+-]\([0-9]*\)/1e-\\1/g"
126: TEST*/