Actual source code: test7.c
slepc-3.12.0 2019-09-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 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: {
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 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 sqrtm(A) - - - - - - -\n");
47: MatView(S,viewer);
48: }
49: /* check error ||S*S-A||_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: MatAXPY(R,-tau,A,SAME_NONZERO_PATTERN);
55: MatNorm(R,NORM_FROBENIUS,&nrm);
56: if (nrm<100*PETSC_MACHINE_EPSILON) {
57: PetscPrintf(PETSC_COMM_WORLD,"||S*S-A||_F < 100*eps\n");
58: } else {
59: PetscPrintf(PETSC_COMM_WORLD,"||S*S-A||_F = %g\n",(double)nrm);
60: }
61: /* check FNEvaluateFunctionMatVec() */
62: MatCreateVecs(A,&v,&f0);
63: MatGetColumnVector(S,f0,0);
64: FNEvaluateFunctionMatVec(fn,A,v);
65: VecAXPY(v,-1.0,f0);
66: VecNorm(v,NORM_2,&nrm);
67: if (nrm>100*PETSC_MACHINE_EPSILON) {
68: PetscPrintf(PETSC_COMM_WORLD,"Warning: the norm of f(A)*e_1-v is %g\n",(double)nrm);
69: }
70: MatDestroy(&S);
71: MatDestroy(&R);
72: VecDestroy(&v);
73: VecDestroy(&f0);
74: return(0);
75: }
77: int main(int argc,char **argv)
78: {
80: FN fn;
81: Mat A;
82: PetscInt i,j,n=10;
83: PetscScalar *As;
84: PetscViewer viewer;
85: PetscBool verbose,inplace;
86: PetscRandom myrand;
87: PetscReal v;
89: SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
90: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
91: PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
92: PetscOptionsHasName(NULL,NULL,"-inplace",&inplace);
93: PetscPrintf(PETSC_COMM_WORLD,"Matrix square root, n=%D.\n",n);
95: /* Create function object */
96: FNCreate(PETSC_COMM_WORLD,&fn);
97: FNSetType(fn,FNSQRT);
98: FNSetFromOptions(fn);
100: /* Set up viewer */
101: PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
102: FNView(fn,viewer);
103: if (verbose) {
104: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);
105: }
107: /* Create matrix */
108: MatCreateSeqDense(PETSC_COMM_SELF,n,n,NULL,&A);
109: PetscObjectSetName((PetscObject)A,"A");
111: /* Compute square root of a symmetric matrix A */
112: MatDenseGetArray(A,&As);
113: for (i=0;i<n;i++) As[i+i*n]=2.5;
114: for (j=1;j<3;j++) {
115: for (i=0;i<n-j;i++) { As[i+(i+j)*n]=1.0; As[(i+j)+i*n]=1.0; }
116: }
117: MatDenseRestoreArray(A,&As);
118: MatSetOption(A,MAT_HERMITIAN,PETSC_TRUE);
119: TestMatSqrt(fn,A,viewer,verbose,inplace);
121: /* Repeat with upper triangular A */
122: MatDenseGetArray(A,&As);
123: for (j=1;j<3;j++) {
124: for (i=0;i<n-j;i++) As[(i+j)+i*n]=0.0;
125: }
126: MatDenseRestoreArray(A,&As);
127: MatSetOption(A,MAT_HERMITIAN,PETSC_FALSE);
128: TestMatSqrt(fn,A,viewer,verbose,inplace);
130: /* Repeat with non-symmetic A */
131: PetscRandomCreate(PETSC_COMM_WORLD,&myrand);
132: PetscRandomSetFromOptions(myrand);
133: PetscRandomSetInterval(myrand,0.0,1.0);
134: MatDenseGetArray(A,&As);
135: for (j=1;j<3;j++) {
136: for (i=0;i<n-j;i++) {
137: PetscRandomGetValueReal(myrand,&v);
138: As[(i+j)+i*n]=v;
139: }
140: }
141: MatDenseRestoreArray(A,&As);
142: PetscRandomDestroy(&myrand);
143: MatSetOption(A,MAT_HERMITIAN,PETSC_FALSE);
144: TestMatSqrt(fn,A,viewer,verbose,inplace);
146: MatDestroy(&A);
147: FNDestroy(&fn);
148: SlepcFinalize();
149: return ierr;
150: }
152: /*TEST
154: test:
155: suffix: 1
156: nsize: 1
157: args: -fn_scale .05,2 -n 100 -fn_method {{0 1 2}shared output}
158: filter: grep -v "computing matrix functions"
159: output_file: output/test7_1.out
160: timeoutfactor: 2
162: test:
163: suffix: 1_sadeghi
164: nsize: 1
165: args: -fn_scale .05,2 -n 100 -fn_method 3
166: requires: !single
167: filter: grep -v "computing matrix functions"
168: output_file: output/test7_1.out
170: test:
171: suffix: 2
172: nsize: 1
173: args: -fn_scale .05,2 -n 100 -inplace -fn_method {{0 1 2}shared output}
174: filter: grep -v "computing matrix functions"
175: output_file: output/test7_1.out
176: timeoutfactor: 2
178: test:
179: suffix: 2_sadeghi
180: nsize: 1
181: args: -fn_scale .05,2 -n 100 -inplace -fn_method 3
182: requires: !single
183: filter: grep -v "computing matrix functions"
184: output_file: output/test7_1.out
186: test:
187: suffix: 3
188: nsize: 3
189: args: -fn_scale .05,2 -n 100 -fn_parallel synchronized
190: filter: grep -v "computing matrix functions" | grep -v "SYNCHRONIZED" | sed -e "s/3 MPI/1 MPI/g"
191: output_file: output/test7_1.out
193: test:
194: suffix: 4
195: nsize: 3
196: args: -fn_scale .05,2 -n 100 -inplace -fn_parallel synchronized
197: filter: grep -v "computing matrix functions" | grep -v "SYNCHRONIZED" | sed -e "s/3 MPI/1 MPI/g"
198: output_file: output/test7_1.out
200: TEST*/