Actual source code: test8.c
slepc-3.16.1 2021-11-17
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2021, 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: FNGetScale(fn,&tau,&eta);
34: /* compute inverse square root */
35: if (inplace) {
36: MatCopy(A,S,SAME_NONZERO_PATTERN);
37: MatIsHermitianKnown(A,&set,&flg);
38: if (set && flg) { MatSetOption(S,MAT_HERMITIAN,PETSC_TRUE); }
39: FNEvaluateFunctionMat(fn,S,NULL);
40: } else {
41: FNEvaluateFunctionMat(fn,A,S);
42: }
43: if (verbose) {
44: PetscPrintf(PETSC_COMM_WORLD,"Matrix A - - - - - - - -\n");
45: MatView(A,viewer);
46: PetscPrintf(PETSC_COMM_WORLD,"Computed inv(sqrtm(A)) - - - - - - -\n");
47: MatView(S,viewer);
48: }
49: /* check error ||S*S*A-I||_F */
50: MatMatMult(S,S,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&R);
51: if (eta!=1.0) {
52: MatScale(R,1.0/(eta*eta));
53: }
54: MatCreateVecs(A,&v,&f0);
55: MatGetColumnVector(S,f0,0);
56: MatCopy(R,S,SAME_NONZERO_PATTERN);
57: MatDestroy(&R);
58: if (tau!=1.0) {
59: MatScale(S,tau);
60: }
61: MatMatMult(S,A,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&R);
62: MatShift(R,-1.0);
63: MatNorm(R,NORM_FROBENIUS,&nrm);
64: MatDestroy(&R);
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: VecDestroy(&v);
79: VecDestroy(&f0);
80: return(0);
81: }
83: int main(int argc,char **argv)
84: {
86: FN fn;
87: Mat A;
88: PetscInt i,j,n=10;
89: PetscScalar x,y,yp,*As;
90: PetscViewer viewer;
91: PetscBool verbose,inplace;
92: PetscRandom myrand;
93: PetscReal v;
94: char strx[50],str[50];
96: SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
97: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
98: PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
99: PetscOptionsHasName(NULL,NULL,"-inplace",&inplace);
100: PetscPrintf(PETSC_COMM_WORLD,"Matrix inverse square root, n=%D.\n",n);
102: /* Create function object */
103: FNCreate(PETSC_COMM_WORLD,&fn);
104: FNSetType(fn,FNINVSQRT);
105: FNSetFromOptions(fn);
107: /* Set up viewer */
108: PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
109: FNView(fn,viewer);
110: if (verbose) {
111: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);
112: }
114: /* Scalar evaluation */
115: x = 2.2;
116: SlepcSNPrintfScalar(strx,sizeof(strx),x,PETSC_FALSE);
117: FNEvaluateFunction(fn,x,&y);
118: FNEvaluateDerivative(fn,x,&yp);
119: SlepcSNPrintfScalar(str,sizeof(str),y,PETSC_FALSE);
120: PetscPrintf(PETSC_COMM_WORLD," f(%s)=%s\n",strx,str);
121: SlepcSNPrintfScalar(str,sizeof(str),yp,PETSC_FALSE);
122: PetscPrintf(PETSC_COMM_WORLD," f'(%s)=%s\n",strx,str);
124: /* Create matrix */
125: MatCreateSeqDense(PETSC_COMM_SELF,n,n,NULL,&A);
126: PetscObjectSetName((PetscObject)A,"A");
128: /* Compute square root of a symmetric matrix A */
129: MatDenseGetArray(A,&As);
130: for (i=0;i<n;i++) As[i+i*n]=2.5;
131: for (j=1;j<3;j++) {
132: for (i=0;i<n-j;i++) { As[i+(i+j)*n]=1.0; As[(i+j)+i*n]=1.0; }
133: }
134: MatDenseRestoreArray(A,&As);
135: MatSetOption(A,MAT_HERMITIAN,PETSC_TRUE);
136: TestMatInvSqrt(fn,A,viewer,verbose,inplace);
138: /* Repeat with upper triangular A */
139: MatDenseGetArray(A,&As);
140: for (j=1;j<3;j++) {
141: for (i=0;i<n-j;i++) As[(i+j)+i*n]=0.0;
142: }
143: MatDenseRestoreArray(A,&As);
144: MatSetOption(A,MAT_HERMITIAN,PETSC_FALSE);
145: TestMatInvSqrt(fn,A,viewer,verbose,inplace);
147: /* Repeat with non-symmetic A */
148: PetscRandomCreate(PETSC_COMM_WORLD,&myrand);
149: PetscRandomSetFromOptions(myrand);
150: PetscRandomSetInterval(myrand,0.0,1.0);
151: MatDenseGetArray(A,&As);
152: for (j=1;j<3;j++) {
153: for (i=0;i<n-j;i++) {
154: PetscRandomGetValueReal(myrand,&v);
155: As[(i+j)+i*n]=v;
156: }
157: }
158: MatDenseRestoreArray(A,&As);
159: PetscRandomDestroy(&myrand);
160: MatSetOption(A,MAT_HERMITIAN,PETSC_FALSE);
161: TestMatInvSqrt(fn,A,viewer,verbose,inplace);
163: MatDestroy(&A);
164: FNDestroy(&fn);
165: SlepcFinalize();
166: return ierr;
167: }
169: /*TEST
171: testset:
172: args: -fn_scale 0.9,0.5 -n 10
173: filter: grep -v "computing matrix functions"
174: output_file: output/test8_1.out
175: test:
176: suffix: 1
177: args: -fn_method {{0 1 2 3}}
178: test:
179: suffix: 1_cuda
180: args: -fn_method 4
181: requires: cuda
182: test:
183: suffix: 1_magma
184: args: -fn_method {{5 6}}
185: requires: cuda magma
186: test:
187: suffix: 2
188: args: -inplace -fn_method {{0 1 2 3}}
189: test:
190: suffix: 2_cuda
191: args: -inplace -fn_method 4
192: requires: cuda
193: test:
194: suffix: 2_magma
195: args: -inplace -fn_method {{5 6}}
196: requires: cuda magma
198: TEST*/