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