Actual source code: test7.c

slepc-3.11.2 2019-07-30
Report Typos and Errors
  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:   MatCreateSeqDense(PETSC_COMM_SELF,n,n,NULL,&R);
 34:   PetscObjectSetName((PetscObject)R,"R");
 35:   FNGetScale(fn,&tau,&eta);
 36:   /* compute 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 sqrtm(A) - - - - - - -\n");
 49:     MatView(S,viewer);
 50:   }
 51:   /* check error ||S*S-A||_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:   MatAXPY(R,-tau,A,SAME_NONZERO_PATTERN);
 57:   MatNorm(R,NORM_FROBENIUS,&nrm);
 58:   if (nrm<100*PETSC_MACHINE_EPSILON) {
 59:     PetscPrintf(PETSC_COMM_WORLD,"||S*S-A||_F < 100*eps\n");
 60:   } else {
 61:     PetscPrintf(PETSC_COMM_WORLD,"||S*S-A||_F = %g\n",(double)nrm);
 62:   }
 63:   /* check FNEvaluateFunctionMatVec() */
 64:   MatCreateVecs(A,&v,&f0);
 65:   MatGetColumnVector(S,f0,0);
 66:   FNEvaluateFunctionMatVec(fn,A,v);
 67:   VecAXPY(v,-1.0,f0);
 68:   VecNorm(v,NORM_2,&nrm);
 69:   if (nrm>100*PETSC_MACHINE_EPSILON) {
 70:     PetscPrintf(PETSC_COMM_WORLD,"Warning: the norm of f(A)*e_1-v is %g\n",(double)nrm);
 71:   }
 72:   MatDestroy(&S);
 73:   MatDestroy(&R);
 74:   VecDestroy(&v);
 75:   VecDestroy(&f0);
 76:   return(0);
 77: }

 79: int main(int argc,char **argv)
 80: {
 82:   FN             fn;
 83:   Mat            A;
 84:   PetscInt       i,j,n=10;
 85:   PetscScalar    *As;
 86:   PetscViewer    viewer;
 87:   PetscBool      verbose,inplace;
 88:   PetscRandom    myrand;
 89:   PetscReal      v;

 91:   SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
 92:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 93:   PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
 94:   PetscOptionsHasName(NULL,NULL,"-inplace",&inplace);
 95:   PetscPrintf(PETSC_COMM_WORLD,"Matrix square root, n=%D.\n",n);

 97:   /* Create function object */
 98:   FNCreate(PETSC_COMM_WORLD,&fn);
 99:   FNSetType(fn,FNSQRT);
100:   FNSetFromOptions(fn);

102:   /* Set up viewer */
103:   PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
104:   FNView(fn,viewer);
105:   if (verbose) {
106:     PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);
107:   }

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:   TestMatSqrt(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:   TestMatSqrt(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:   TestMatSqrt(fn,A,viewer,verbose,inplace);

148:   MatDestroy(&A);
149:   FNDestroy(&fn);
150:   SlepcFinalize();
151:   return ierr;
152: }

154: /*TEST

156:    test:
157:       suffix: 1
158:       nsize: 1
159:       args: -fn_scale .05,2 -n 100 -fn_method {{0 1 2}shared output}
160:       filter: grep -v "computing matrix functions"
161:       output_file: output/test7_1.out
162:       timeoutfactor: 2

164:    test:
165:       suffix: 1_sadeghi
166:       nsize: 1
167:       args: -fn_scale .05,2 -n 100 -fn_method 3
168:       requires: !single
169:       filter: grep -v "computing matrix functions"
170:       output_file: output/test7_1.out

172:    test:
173:       suffix: 2
174:       nsize: 1
175:       args: -fn_scale .05,2 -n 100 -inplace -fn_method {{0 1 2}shared output}
176:       filter: grep -v "computing matrix functions"
177:       output_file: output/test7_1.out
178:       timeoutfactor: 2

180:    test:
181:       suffix: 2_sadeghi
182:       nsize: 1
183:       args: -fn_scale .05,2 -n 100 -inplace -fn_method 3
184:       requires: !single
185:       filter: grep -v "computing matrix functions"
186:       output_file: output/test7_1.out

188:    test:
189:       suffix: 3
190:       nsize: 3
191:       args: -fn_scale .05,2 -n 100 -fn_parallel synchronized
192:       filter: grep -v "computing matrix functions" | grep -v "SYNCHRONIZED" | sed -e "s/3 MPI/1 MPI/g"
193:       output_file: output/test7_1.out

195:    test:
196:       suffix: 4
197:       nsize: 3
198:       args: -fn_scale .05,2 -n 100 -inplace -fn_parallel synchronized
199:       filter: grep -v "computing matrix functions" | grep -v "SYNCHRONIZED" | sed -e "s/3 MPI/1 MPI/g"
200:       output_file: output/test7_1.out

202: TEST*/