Actual source code: test8.c
slepc-3.11.2 2019-07-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 matrix inverse square root.\n\n";
13: #include <slepcfn.h>
15: /*
16: Compute matrix inverse square root B = inv(sqrtm(A))
17: Check result as norm(B*B*A-I)
18: */
19: PetscErrorCode TestMatInvSqrt(FN fn,Mat A,PetscViewer viewer,PetscBool verbose,PetscBool inplace)
20: {
22: PetscScalar tau,eta;
23: PetscReal nrm;
24: PetscBool set,flg;
25: PetscInt n;
26: Mat S,R;
27: Vec v,f0;
30: MatGetSize(A,&n,NULL);
31: MatCreateSeqDense(PETSC_COMM_SELF,n,n,NULL,&S);
32: PetscObjectSetName((PetscObject)S,"S");
33: MatCreateSeqDense(PETSC_COMM_SELF,n,n,NULL,&R);
34: PetscObjectSetName((PetscObject)R,"R");
35: FNGetScale(fn,&tau,&eta);
36: /* compute inverse square root */
37: if (inplace) {
38: MatCopy(A,S,SAME_NONZERO_PATTERN);
39: MatIsHermitianKnown(A,&set,&flg);
40: if (set && flg) { MatSetOption(S,MAT_HERMITIAN,PETSC_TRUE); }
41: FNEvaluateFunctionMat(fn,S,NULL);
42: } else {
43: FNEvaluateFunctionMat(fn,A,S);
44: }
45: if (verbose) {
46: PetscPrintf(PETSC_COMM_WORLD,"Matrix A - - - - - - - -\n");
47: MatView(A,viewer);
48: PetscPrintf(PETSC_COMM_WORLD,"Computed inv(sqrtm(A)) - - - - - - -\n");
49: MatView(S,viewer);
50: }
51: /* check error ||S*S*A-I||_F */
52: MatMatMult(S,S,MAT_REUSE_MATRIX,PETSC_DEFAULT,&R);
53: if (eta!=1.0) {
54: MatScale(R,1.0/(eta*eta));
55: }
56: MatCreateVecs(A,&v,&f0);
57: MatGetColumnVector(S,f0,0);
58: MatCopy(R,S,SAME_NONZERO_PATTERN);
59: if (tau!=1.0) {
60: MatScale(S,tau);
61: }
62: MatMatMult(S,A,MAT_REUSE_MATRIX,PETSC_DEFAULT,&R);
63: MatShift(R,-1.0);
64: MatNorm(R,NORM_FROBENIUS,&nrm);
65: if (nrm<100*PETSC_MACHINE_EPSILON) {
66: PetscPrintf(PETSC_COMM_WORLD,"||S*S*A-I||_F < 100*eps\n");
67: } else {
68: PetscPrintf(PETSC_COMM_WORLD,"||S*S*A-I||_F = %g\n",(double)nrm);
69: }
70: /* check FNEvaluateFunctionMatVec() */
71: FNEvaluateFunctionMatVec(fn,A,v);
72: VecAXPY(v,-1.0,f0);
73: VecNorm(v,NORM_2,&nrm);
74: if (nrm>100*PETSC_MACHINE_EPSILON) {
75: PetscPrintf(PETSC_COMM_WORLD,"Warning: the norm of f(A)*e_1-v is %g\n",(double)nrm);
76: }
77: MatDestroy(&S);
78: MatDestroy(&R);
79: VecDestroy(&v);
80: VecDestroy(&f0);
81: return(0);
82: }
84: int main(int argc,char **argv)
85: {
87: FN fn;
88: Mat A;
89: PetscInt i,j,n=10;
90: PetscScalar x,y,yp,*As;
91: PetscViewer viewer;
92: PetscBool verbose,inplace;
93: PetscRandom myrand;
94: PetscReal v;
95: char strx[50],str[50];
97: SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
98: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
99: PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
100: PetscOptionsHasName(NULL,NULL,"-inplace",&inplace);
101: PetscPrintf(PETSC_COMM_WORLD,"Matrix inverse square root, n=%D.\n",n);
103: /* Create function object */
104: FNCreate(PETSC_COMM_WORLD,&fn);
105: FNSetType(fn,FNINVSQRT);
106: FNSetFromOptions(fn);
108: /* Set up viewer */
109: PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
110: FNView(fn,viewer);
111: if (verbose) {
112: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);
113: }
115: /* Scalar evaluation */
116: x = 2.2;
117: SlepcSNPrintfScalar(strx,50,x,PETSC_FALSE);
118: FNEvaluateFunction(fn,x,&y);
119: FNEvaluateDerivative(fn,x,&yp);
120: SlepcSNPrintfScalar(str,50,y,PETSC_FALSE);
121: PetscPrintf(PETSC_COMM_WORLD," f(%s)=%s\n",strx,str);
122: SlepcSNPrintfScalar(str,50,yp,PETSC_FALSE);
123: PetscPrintf(PETSC_COMM_WORLD," f'(%s)=%s\n",strx,str);
125: /* Create matrix */
126: MatCreateSeqDense(PETSC_COMM_SELF,n,n,NULL,&A);
127: PetscObjectSetName((PetscObject)A,"A");
129: /* Compute square root of a symmetric matrix A */
130: MatDenseGetArray(A,&As);
131: for (i=0;i<n;i++) As[i+i*n]=2.5;
132: for (j=1;j<3;j++) {
133: for (i=0;i<n-j;i++) { As[i+(i+j)*n]=1.0; As[(i+j)+i*n]=1.0; }
134: }
135: MatDenseRestoreArray(A,&As);
136: MatSetOption(A,MAT_HERMITIAN,PETSC_TRUE);
137: TestMatInvSqrt(fn,A,viewer,verbose,inplace);
139: /* Repeat with upper triangular A */
140: MatDenseGetArray(A,&As);
141: for (j=1;j<3;j++) {
142: for (i=0;i<n-j;i++) As[(i+j)+i*n]=0.0;
143: }
144: MatDenseRestoreArray(A,&As);
145: MatSetOption(A,MAT_HERMITIAN,PETSC_FALSE);
146: TestMatInvSqrt(fn,A,viewer,verbose,inplace);
148: /* Repeat with non-symmetic A */
149: PetscRandomCreate(PETSC_COMM_WORLD,&myrand);
150: PetscRandomSetFromOptions(myrand);
151: PetscRandomSetInterval(myrand,0.0,1.0);
152: MatDenseGetArray(A,&As);
153: for (j=1;j<3;j++) {
154: for (i=0;i<n-j;i++) {
155: PetscRandomGetValueReal(myrand,&v);
156: As[(i+j)+i*n]=v;
157: }
158: }
159: MatDenseRestoreArray(A,&As);
160: PetscRandomDestroy(&myrand);
161: MatSetOption(A,MAT_HERMITIAN,PETSC_FALSE);
162: TestMatInvSqrt(fn,A,viewer,verbose,inplace);
164: MatDestroy(&A);
165: FNDestroy(&fn);
166: SlepcFinalize();
167: return ierr;
168: }
170: /*TEST
172: test:
173: suffix: 1
174: nsize: 1
175: args: -fn_scale 0.9,0.5 -n 10 -fn_method {{0 1 2}shared output}
176: filter: grep -v "computing matrix functions"
177: output_file: output/test8_1.out
179: test:
180: suffix: 2
181: nsize: 1
182: args: -fn_scale 0.9,0.5 -n 10 -inplace -fn_method {{0 1 2}shared output}
183: filter: grep -v "computing matrix functions"
184: output_file: output/test8_1.out
186: TEST*/