Actual source code: test1.c
slepc-3.22.1 2024-10-28
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 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";
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: PetscBool flg,testprefix=PETSC_FALSE,viewmatrices=PETSC_FALSE;
26: const char *prefix;
27: LMEType type;
28: LMEProblemType ptype;
29: PetscViewerAndFormat *vf;
31: PetscFunctionBeginUser;
32: PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
34: PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
35: PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,&flg));
36: if (!flg) m=n;
37: N = n*m;
38: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nLyapunov equation, N=%" PetscInt_FMT " (%" PetscInt_FMT "x%" PetscInt_FMT " grid)\n\n",N,n,m));
39: PetscCall(PetscOptionsGetBool(NULL,NULL,"-test_prefix",&testprefix,NULL));
40: PetscCall(PetscOptionsGetBool(NULL,NULL,"-view_matrices",&viewmatrices,NULL));
42: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
43: Create the 2-D Laplacian, A
44: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
46: PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
47: PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N));
48: PetscCall(MatSetFromOptions(A));
49: PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
50: for (II=Istart;II<Iend;II++) {
51: i = II/n; j = II-i*n;
52: if (i>0) PetscCall(MatSetValue(A,II,II-n,1.0,INSERT_VALUES));
53: if (i<m-1) PetscCall(MatSetValue(A,II,II+n,1.0,INSERT_VALUES));
54: if (j>0) PetscCall(MatSetValue(A,II,II-1,1.0,INSERT_VALUES));
55: if (j<n-1) PetscCall(MatSetValue(A,II,II+1,1.0,INSERT_VALUES));
56: PetscCall(MatSetValue(A,II,II,-4.0,INSERT_VALUES));
57: }
58: PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
59: PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
61: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
62: Create a low-rank Mat to store the right-hand side C = C1*C1'
63: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
65: PetscCall(MatCreate(PETSC_COMM_WORLD,&C1));
66: PetscCall(MatSetSizes(C1,PETSC_DECIDE,PETSC_DECIDE,N,2));
67: PetscCall(MatSetType(C1,MATDENSE));
68: PetscCall(MatGetOwnershipRange(C1,&Istart,&Iend));
69: PetscCall(MatDenseGetArray(C1,&u));
70: for (i=Istart;i<Iend;i++) {
71: if (i<N/2) u[i-Istart] = 1.0;
72: if (i==0) u[i+Iend-2*Istart] = -2.0;
73: if (i==1) u[i+Iend-2*Istart] = -1.0;
74: if (i==2) u[i+Iend-2*Istart] = -1.0;
75: }
76: PetscCall(MatDenseRestoreArray(C1,&u));
77: PetscCall(MatAssemblyBegin(C1,MAT_FINAL_ASSEMBLY));
78: PetscCall(MatAssemblyEnd(C1,MAT_FINAL_ASSEMBLY));
79: PetscCall(MatCreateLRC(NULL,C1,NULL,NULL,&C));
80: PetscCall(MatDestroy(&C1));
82: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
83: Create the solver and set various options
84: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
85: PetscCall(LMECreate(PETSC_COMM_WORLD,&lme));
86: PetscCall(LMESetProblemType(lme,LME_SYLVESTER));
87: PetscCall(LMEGetProblemType(lme,&ptype));
88: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Equation type set to %d\n",ptype));
89: PetscCall(LMESetProblemType(lme,LME_LYAPUNOV));
90: PetscCall(LMEGetProblemType(lme,&ptype));
91: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Equation type changed to %d\n",ptype));
92: PetscCall(LMESetCoefficients(lme,A,NULL,NULL,NULL));
93: PetscCall(LMESetRHS(lme,C));
95: /* test prefix usage */
96: if (testprefix) {
97: PetscCall(LMESetOptionsPrefix(lme,"check_"));
98: PetscCall(LMEAppendOptionsPrefix(lme,"myprefix_"));
99: PetscCall(LMEGetOptionsPrefix(lme,&prefix));
100: PetscCall(PetscPrintf(PETSC_COMM_WORLD," LME prefix is currently: %s\n",prefix));
101: }
103: /* test some interface functions */
104: PetscCall(LMEGetCoefficients(lme,&B,NULL,NULL,NULL));
105: if (viewmatrices) PetscCall(MatView(B,PETSC_VIEWER_STDOUT_WORLD));
106: PetscCall(LMEGetRHS(lme,&D));
107: if (viewmatrices) PetscCall(MatView(D,PETSC_VIEWER_STDOUT_WORLD));
108: PetscCall(LMESetTolerances(lme,PETSC_CURRENT,100));
109: PetscCall(LMESetDimensions(lme,21));
110: PetscCall(LMESetErrorIfNotConverged(lme,PETSC_TRUE));
111: /* test monitors */
112: PetscCall(PetscViewerAndFormatCreate(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_DEFAULT,&vf));
113: PetscCall(LMEMonitorSet(lme,(PetscErrorCode (*)(LME,PetscInt,PetscReal,void*))LMEMonitorDefault,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy));
114: /* PetscCall(LMEMonitorCancel(lme)); */
115: PetscCall(LMESetFromOptions(lme));
117: PetscCall(LMEGetType(lme,&type));
118: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Solver being used: %s\n",type));
120: /* query properties and print them */
121: PetscCall(LMEGetTolerances(lme,&tol,&maxit));
122: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Tolerance: %g, max iterations: %" PetscInt_FMT "\n",(double)tol,maxit));
123: PetscCall(LMEGetDimensions(lme,&ncv));
124: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Subspace dimension: %" PetscInt_FMT "\n",ncv));
125: PetscCall(LMEGetErrorIfNotConverged(lme,&flg));
126: if (flg) PetscCall(PetscPrintf(PETSC_COMM_WORLD," Erroring out if convergence fails\n"));
128: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
129: Solve the matrix equation and compute residual error
130: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
132: PetscCall(LMESolve(lme));
133: PetscCall(LMEGetErrorEstimate(lme,&errest));
134: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Error estimate reported by the solver: %.4g\n",(double)errest));
135: PetscCall(LMEComputeError(lme,&error));
136: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Computed residual norm: %.4g\n\n",(double)error));
138: /*
139: Free work space
140: */
141: PetscCall(LMEDestroy(&lme));
142: PetscCall(MatDestroy(&A));
143: PetscCall(MatDestroy(&C));
144: PetscCall(SlepcFinalize());
145: return 0;
146: }
148: /*TEST
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/"
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/"
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
168: TEST*/