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