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 LME interface functions, based on ex32.c.\n\n"
12 : "The command line options are:\n"
13 : " -n <n>, where <n> = number of grid subdivisions in x dimension.\n"
14 : " -m <m>, where <m> = number of grid subdivisions in y dimension.\n\n";
15 :
16 : #include <slepclme.h>
17 :
18 2 : int main(int argc,char **argv)
19 : {
20 2 : Mat A,B,C,C1,D;
21 2 : LME lme;
22 2 : PetscReal tol,errest,error;
23 2 : PetscScalar *u;
24 2 : PetscInt N,n=10,m,Istart,Iend,II,maxit,ncv,i,j;
25 2 : PetscBool flg,testprefix=PETSC_FALSE,viewmatrices=PETSC_FALSE;
26 2 : const char *prefix;
27 2 : LMEType type;
28 2 : LMEProblemType ptype;
29 2 : PetscViewerAndFormat *vf;
30 :
31 2 : PetscFunctionBeginUser;
32 2 : PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
33 :
34 2 : PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
35 2 : PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,&flg));
36 2 : if (!flg) m=n;
37 2 : N = n*m;
38 2 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nLyapunov equation, N=%" PetscInt_FMT " (%" PetscInt_FMT "x%" PetscInt_FMT " grid)\n\n",N,n,m));
39 2 : PetscCall(PetscOptionsGetBool(NULL,NULL,"-test_prefix",&testprefix,NULL));
40 2 : PetscCall(PetscOptionsGetBool(NULL,NULL,"-view_matrices",&viewmatrices,NULL));
41 :
42 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
43 : Create the 2-D Laplacian, A
44 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
45 :
46 2 : PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
47 2 : PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N));
48 2 : PetscCall(MatSetFromOptions(A));
49 2 : PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
50 202 : for (II=Istart;II<Iend;II++) {
51 200 : i = II/n; j = II-i*n;
52 200 : if (i>0) PetscCall(MatSetValue(A,II,II-n,1.0,INSERT_VALUES));
53 200 : if (i<m-1) PetscCall(MatSetValue(A,II,II+n,1.0,INSERT_VALUES));
54 200 : if (j>0) PetscCall(MatSetValue(A,II,II-1,1.0,INSERT_VALUES));
55 200 : if (j<n-1) PetscCall(MatSetValue(A,II,II+1,1.0,INSERT_VALUES));
56 200 : PetscCall(MatSetValue(A,II,II,-4.0,INSERT_VALUES));
57 : }
58 2 : PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
59 2 : PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
60 :
61 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
62 : Create a low-rank Mat to store the right-hand side C = C1*C1'
63 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
64 :
65 2 : PetscCall(MatCreate(PETSC_COMM_WORLD,&C1));
66 2 : PetscCall(MatSetSizes(C1,PETSC_DECIDE,PETSC_DECIDE,N,2));
67 2 : PetscCall(MatSetType(C1,MATDENSE));
68 2 : PetscCall(MatGetOwnershipRange(C1,&Istart,&Iend));
69 2 : PetscCall(MatDenseGetArray(C1,&u));
70 202 : for (i=Istart;i<Iend;i++) {
71 200 : if (i<N/2) u[i-Istart] = 1.0;
72 200 : if (i==0) u[i+Iend-2*Istart] = -2.0;
73 200 : if (i==1) u[i+Iend-2*Istart] = -1.0;
74 200 : if (i==2) u[i+Iend-2*Istart] = -1.0;
75 : }
76 2 : PetscCall(MatDenseRestoreArray(C1,&u));
77 2 : PetscCall(MatAssemblyBegin(C1,MAT_FINAL_ASSEMBLY));
78 2 : PetscCall(MatAssemblyEnd(C1,MAT_FINAL_ASSEMBLY));
79 2 : PetscCall(MatCreateLRC(NULL,C1,NULL,NULL,&C));
80 2 : PetscCall(MatDestroy(&C1));
81 :
82 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
83 : Create the solver and set various options
84 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
85 2 : PetscCall(LMECreate(PETSC_COMM_WORLD,&lme));
86 2 : PetscCall(LMESetProblemType(lme,LME_SYLVESTER));
87 2 : PetscCall(LMEGetProblemType(lme,&ptype));
88 2 : PetscCall(PetscPrintf(PETSC_COMM_WORLD," Equation type set to %d\n",ptype));
89 2 : PetscCall(LMESetProblemType(lme,LME_LYAPUNOV));
90 2 : PetscCall(LMEGetProblemType(lme,&ptype));
91 2 : PetscCall(PetscPrintf(PETSC_COMM_WORLD," Equation type changed to %d\n",ptype));
92 2 : PetscCall(LMESetCoefficients(lme,A,NULL,NULL,NULL));
93 2 : PetscCall(LMESetRHS(lme,C));
94 :
95 : /* test prefix usage */
96 2 : if (testprefix) {
97 1 : PetscCall(LMESetOptionsPrefix(lme,"check_"));
98 1 : PetscCall(LMEAppendOptionsPrefix(lme,"myprefix_"));
99 1 : PetscCall(LMEGetOptionsPrefix(lme,&prefix));
100 1 : PetscCall(PetscPrintf(PETSC_COMM_WORLD," LME prefix is currently: %s\n",prefix));
101 : }
102 :
103 : /* test some interface functions */
104 2 : PetscCall(LMEGetCoefficients(lme,&B,NULL,NULL,NULL));
105 2 : if (viewmatrices) PetscCall(MatView(B,PETSC_VIEWER_STDOUT_WORLD));
106 2 : PetscCall(LMEGetRHS(lme,&D));
107 2 : if (viewmatrices) PetscCall(MatView(D,PETSC_VIEWER_STDOUT_WORLD));
108 2 : PetscCall(LMESetTolerances(lme,PETSC_CURRENT,100));
109 2 : PetscCall(LMESetDimensions(lme,21));
110 2 : PetscCall(LMESetErrorIfNotConverged(lme,PETSC_TRUE));
111 : /* test monitors */
112 2 : PetscCall(PetscViewerAndFormatCreate(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_DEFAULT,&vf));
113 2 : PetscCall(LMEMonitorSet(lme,(PetscErrorCode (*)(LME,PetscInt,PetscReal,void*))LMEMonitorDefault,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy));
114 : /* PetscCall(LMEMonitorCancel(lme)); */
115 2 : PetscCall(LMESetFromOptions(lme));
116 :
117 2 : PetscCall(LMEGetType(lme,&type));
118 2 : PetscCall(PetscPrintf(PETSC_COMM_WORLD," Solver being used: %s\n",type));
119 :
120 : /* query properties and print them */
121 2 : PetscCall(LMEGetTolerances(lme,&tol,&maxit));
122 2 : PetscCall(PetscPrintf(PETSC_COMM_WORLD," Tolerance: %g, max iterations: %" PetscInt_FMT "\n",(double)tol,maxit));
123 2 : PetscCall(LMEGetDimensions(lme,&ncv));
124 2 : PetscCall(PetscPrintf(PETSC_COMM_WORLD," Subspace dimension: %" PetscInt_FMT "\n",ncv));
125 2 : PetscCall(LMEGetErrorIfNotConverged(lme,&flg));
126 2 : if (flg) PetscCall(PetscPrintf(PETSC_COMM_WORLD," Erroring out if convergence fails\n"));
127 :
128 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
129 : Solve the matrix equation and compute residual error
130 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
131 :
132 2 : PetscCall(LMESolve(lme));
133 2 : PetscCall(LMEGetErrorEstimate(lme,&errest));
134 2 : PetscCall(PetscPrintf(PETSC_COMM_WORLD," Error estimate reported by the solver: %.4g\n",(double)errest));
135 2 : PetscCall(LMEComputeError(lme,&error));
136 2 : PetscCall(PetscPrintf(PETSC_COMM_WORLD," Computed residual norm: %.4g\n\n",(double)error));
137 :
138 : /*
139 : Free work space
140 : */
141 2 : PetscCall(LMEDestroy(&lme));
142 2 : PetscCall(MatDestroy(&A));
143 2 : PetscCall(MatDestroy(&C));
144 2 : PetscCall(SlepcFinalize());
145 : return 0;
146 : }
147 :
148 : /*TEST
149 :
150 : test:
151 : suffix: 1
152 : args: -lme_monitor_cancel -lme_converged_reason -lme_view -view_matrices -log_exclude lme,bv
153 : requires: double
154 : filter: sed -e "s/4.0[0-9]*e-10/4.03e-10/"
155 :
156 : test:
157 : suffix: 2
158 : args: -test_prefix -check_myprefix_lme_monitor
159 : requires: double
160 : filter: sed -e "s/estimate [0-9]\.[0-9]*e[+-]\([0-9]*\)/estimate (removed)/g" | sed -e "s/4.0[0-9]*e-10/4.03e-10/"
161 :
162 : test:
163 : suffix: 3
164 : args: -lme_monitor_cancel -info -lme_monitor draw::draw_lg -draw_virtual
165 : requires: x double
166 : filter: sed -e "s/equation = [0-9]\.[0-9]*e[+-]\([0-9]*\)/equation = (removed)/g" | sed -e "s/4.0[0-9]*e-10/4.03e-10/" | grep -v Comm | grep -v machine | grep -v PetscGetHostName | grep -v OpenMP | grep -v Colormap | grep -v "Rank of the Cholesky factor" | grep -v "potrf failed" | grep -v "querying" | grep -v FPTrap | grep -v Device | grep -v SetNumThread
167 :
168 : TEST*/
|