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