Actual source code: test7.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 square root.\n\n";
13: #include <slepcfn.h>
15: /*
16: Compute matrix square root B = sqrtm(A)
17: Check result as norm(B*B-A)
18: */
19: PetscErrorCode TestMatSqrt(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 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 sqrtm(A) - - - - - - -\n");
44: MatView(S,viewer);
45: }
46: /* check error ||S*S-A||_F */
47: MatMatMult(S,S,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&R);
48: if (eta!=1.0) MatScale(R,1.0/(eta*eta));
49: MatAXPY(R,-tau,A,SAME_NONZERO_PATTERN);
50: MatNorm(R,NORM_FROBENIUS,&nrm);
51: if (nrm<100*PETSC_MACHINE_EPSILON) PetscPrintf(PETSC_COMM_WORLD,"||S*S-A||_F < 100*eps\n");
52: else PetscPrintf(PETSC_COMM_WORLD,"||S*S-A||_F = %g\n",(double)nrm);
53: /* check FNEvaluateFunctionMatVec() */
54: MatCreateVecs(A,&v,&f0);
55: MatGetColumnVector(S,f0,0);
56: FNEvaluateFunctionMatVec(fn,A,v);
57: VecAXPY(v,-1.0,f0);
58: VecNorm(v,NORM_2,&nrm);
59: if (nrm>100*PETSC_MACHINE_EPSILON) PetscPrintf(PETSC_COMM_WORLD,"Warning: the norm of f(A)*e_1-v is %g\n",(double)nrm);
60: MatDestroy(&S);
61: MatDestroy(&R);
62: VecDestroy(&v);
63: VecDestroy(&f0);
64: PetscFunctionReturn(0);
65: }
67: int main(int argc,char **argv)
68: {
69: FN fn;
70: Mat A;
71: PetscInt i,j,n=10;
72: PetscScalar *As;
73: PetscViewer viewer;
74: PetscBool verbose,inplace;
75: PetscRandom myrand;
76: PetscReal v;
78: SlepcInitialize(&argc,&argv,(char*)0,help);
79: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
80: PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
81: PetscOptionsHasName(NULL,NULL,"-inplace",&inplace);
82: PetscPrintf(PETSC_COMM_WORLD,"Matrix square root, n=%" PetscInt_FMT ".\n",n);
84: /* Create function object */
85: FNCreate(PETSC_COMM_WORLD,&fn);
86: FNSetType(fn,FNSQRT);
87: FNSetFromOptions(fn);
89: /* Set up viewer */
90: PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
91: FNView(fn,viewer);
92: if (verbose) PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);
94: /* Create matrix */
95: MatCreateSeqDense(PETSC_COMM_SELF,n,n,NULL,&A);
96: PetscObjectSetName((PetscObject)A,"A");
98: /* Compute square root of a symmetric matrix A */
99: MatDenseGetArray(A,&As);
100: for (i=0;i<n;i++) As[i+i*n]=2.5;
101: for (j=1;j<3;j++) {
102: for (i=0;i<n-j;i++) { As[i+(i+j)*n]=1.0; As[(i+j)+i*n]=1.0; }
103: }
104: MatDenseRestoreArray(A,&As);
105: MatSetOption(A,MAT_HERMITIAN,PETSC_TRUE);
106: TestMatSqrt(fn,A,viewer,verbose,inplace);
108: /* Repeat with upper triangular A */
109: MatDenseGetArray(A,&As);
110: for (j=1;j<3;j++) {
111: for (i=0;i<n-j;i++) As[(i+j)+i*n]=0.0;
112: }
113: MatDenseRestoreArray(A,&As);
114: MatSetOption(A,MAT_HERMITIAN,PETSC_FALSE);
115: TestMatSqrt(fn,A,viewer,verbose,inplace);
117: /* Repeat with non-symmetic A */
118: PetscRandomCreate(PETSC_COMM_WORLD,&myrand);
119: PetscRandomSetFromOptions(myrand);
120: PetscRandomSetInterval(myrand,0.0,1.0);
121: MatDenseGetArray(A,&As);
122: for (j=1;j<3;j++) {
123: for (i=0;i<n-j;i++) {
124: PetscRandomGetValueReal(myrand,&v);
125: As[(i+j)+i*n]=v;
126: }
127: }
128: MatDenseRestoreArray(A,&As);
129: PetscRandomDestroy(&myrand);
130: MatSetOption(A,MAT_HERMITIAN,PETSC_FALSE);
131: TestMatSqrt(fn,A,viewer,verbose,inplace);
133: MatDestroy(&A);
134: FNDestroy(&fn);
135: SlepcFinalize();
136: return 0;
137: }
139: /*TEST
141: testset:
142: args: -fn_scale .05,2 -n 100
143: filter: grep -v "computing matrix functions"
144: output_file: output/test7_1.out
145: timeoutfactor: 2
146: test:
147: suffix: 1
148: args: -fn_method {{0 1 2}}
149: test:
150: suffix: 1_sadeghi
151: args: -fn_method 3
152: requires: !single
153: test:
154: suffix: 1_cuda
155: args: -fn_method 4
156: requires: cuda !single
157: test:
158: suffix: 1_magma
159: args: -fn_method {{5 6}}
160: requires: cuda magma !single
161: test:
162: suffix: 2
163: args: -inplace -fn_method {{0 1 2}}
164: test:
165: suffix: 2_sadeghi
166: args: -inplace -fn_method 3
167: requires: !single
168: test:
169: suffix: 2_cuda
170: args: -inplace -fn_method 4
171: requires: cuda !single
172: test:
173: suffix: 2_magma
174: args: -inplace -fn_method {{5 6}}
175: requires: cuda magma !single
177: testset:
178: nsize: 3
179: args: -fn_scale .05,2 -n 100 -fn_parallel synchronized
180: filter: grep -v "computing matrix functions" | grep -v "SYNCHRONIZED" | sed -e "s/3 MPI/1 MPI/g"
181: output_file: output/test7_1.out
182: test:
183: suffix: 3
184: test:
185: suffix: 3_inplace
186: args: -inplace
188: TEST*/