Actual source code: test8.c

slepc-3.18.2 2023-01-26
Report Typos and Errors
  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,Acopy;
 26:   Vec            v,f0;

 29:   MatGetSize(A,&n,NULL);
 30:   MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&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 {
 40:     MatDuplicate(A,MAT_COPY_VALUES,&Acopy);
 41:     FNEvaluateFunctionMat(fn,A,S);
 42:     /* check that A has not been modified */
 43:     MatAXPY(Acopy,-1.0,A,SAME_NONZERO_PATTERN);
 44:     MatNorm(Acopy,NORM_1,&nrm);
 45:     if (nrm>100*PETSC_MACHINE_EPSILON) PetscPrintf(PETSC_COMM_WORLD,"Warning: the input matrix has changed by %g\n",(double)nrm);
 46:     MatDestroy(&Acopy);
 47:   }
 48:   if (verbose) {
 49:     PetscPrintf(PETSC_COMM_WORLD,"Matrix A - - - - - - - -\n");
 50:     MatView(A,viewer);
 51:     PetscPrintf(PETSC_COMM_WORLD,"Computed inv(sqrtm(A)) - - - - - - -\n");
 52:     MatView(S,viewer);
 53:   }
 54:   /* check error ||S*S*A-I||_F */
 55:   MatMatMult(S,S,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&R);
 56:   if (eta!=1.0) MatScale(R,1.0/(eta*eta));
 57:   MatCreateVecs(A,&v,&f0);
 58:   MatGetColumnVector(S,f0,0);
 59:   MatCopy(R,S,SAME_NONZERO_PATTERN);
 60:   MatDestroy(&R);
 61:   if (tau!=1.0) MatScale(S,tau);
 62:   MatMatMult(S,A,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&R);
 63:   MatShift(R,-1.0);
 64:   MatNorm(R,NORM_FROBENIUS,&nrm);
 65:   MatDestroy(&R);
 66:   if (nrm<100*PETSC_MACHINE_EPSILON) PetscPrintf(PETSC_COMM_WORLD,"||S*S*A-I||_F < 100*eps\n");
 67:   else PetscPrintf(PETSC_COMM_WORLD,"||S*S*A-I||_F = %g\n",(double)nrm);
 68:   /* check FNEvaluateFunctionMatVec() */
 69:   FNEvaluateFunctionMatVec(fn,A,v);
 70:   VecAXPY(v,-1.0,f0);
 71:   VecNorm(v,NORM_2,&nrm);
 72:   if (nrm>100*PETSC_MACHINE_EPSILON) PetscPrintf(PETSC_COMM_WORLD,"Warning: the norm of f(A)*e_1-v is %g\n",(double)nrm);
 73:   MatDestroy(&S);
 74:   VecDestroy(&v);
 75:   VecDestroy(&f0);
 76:   return 0;
 77: }

 79: int main(int argc,char **argv)
 80: {
 81:   FN             fn;
 82:   Mat            A=NULL;
 83:   PetscInt       i,j,n=10;
 84:   PetscScalar    x,y,yp,*As;
 85:   PetscViewer    viewer;
 86:   PetscBool      verbose,inplace,matcuda;
 87:   PetscRandom    myrand;
 88:   PetscReal      v;
 89:   char           strx[50],str[50];

 92:   SlepcInitialize(&argc,&argv,(char*)0,help);
 93:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 94:   PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
 95:   PetscOptionsHasName(NULL,NULL,"-inplace",&inplace);
 96:   PetscOptionsHasName(NULL,NULL,"-matcuda",&matcuda);
 97:   PetscPrintf(PETSC_COMM_WORLD,"Matrix inverse square root, n=%" PetscInt_FMT ".\n",n);

 99:   /* Create function object */
100:   FNCreate(PETSC_COMM_WORLD,&fn);
101:   FNSetType(fn,FNINVSQRT);
102:   FNSetFromOptions(fn);

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

109:   /* Scalar evaluation */
110:   x = 2.2;
111:   SlepcSNPrintfScalar(strx,sizeof(strx),x,PETSC_FALSE);
112:   FNEvaluateFunction(fn,x,&y);
113:   FNEvaluateDerivative(fn,x,&yp);
114:   SlepcSNPrintfScalar(str,sizeof(str),y,PETSC_FALSE);
115:   PetscPrintf(PETSC_COMM_WORLD,"  f(%s)=%s\n",strx,str);
116:   SlepcSNPrintfScalar(str,sizeof(str),yp,PETSC_FALSE);
117:   PetscPrintf(PETSC_COMM_WORLD,"  f'(%s)=%s\n",strx,str);

119:   /* Create matrix */
120:   if (matcuda) {
121: #if defined(PETSC_HAVE_CUDA)
122:     MatCreateSeqDenseCUDA(PETSC_COMM_SELF,n,n,NULL,&A);
123: #endif
124:   } else MatCreateSeqDense(PETSC_COMM_SELF,n,n,NULL,&A);
125:   PetscObjectSetName((PetscObject)A,"A");

127:   /* Compute square root of a symmetric matrix A */
128:   MatDenseGetArray(A,&As);
129:   for (i=0;i<n;i++) As[i+i*n]=2.5;
130:   for (j=1;j<3;j++) {
131:     for (i=0;i<n-j;i++) { As[i+(i+j)*n]=1.0; As[(i+j)+i*n]=1.0; }
132:   }
133:   MatDenseRestoreArray(A,&As);
134:   MatSetOption(A,MAT_HERMITIAN,PETSC_TRUE);
135:   TestMatInvSqrt(fn,A,viewer,verbose,inplace);

137:   /* Repeat with upper triangular A */
138:   MatDenseGetArray(A,&As);
139:   for (j=1;j<3;j++) {
140:     for (i=0;i<n-j;i++) As[(i+j)+i*n]=0.0;
141:   }
142:   MatDenseRestoreArray(A,&As);
143:   MatSetOption(A,MAT_HERMITIAN,PETSC_FALSE);
144:   TestMatInvSqrt(fn,A,viewer,verbose,inplace);

146:   /* Repeat with non-symmetic A */
147:   PetscRandomCreate(PETSC_COMM_WORLD,&myrand);
148:   PetscRandomSetFromOptions(myrand);
149:   PetscRandomSetInterval(myrand,0.0,1.0);
150:   MatDenseGetArray(A,&As);
151:   for (j=1;j<3;j++) {
152:     for (i=0;i<n-j;i++) {
153:       PetscRandomGetValueReal(myrand,&v);
154:       As[(i+j)+i*n]=v;
155:     }
156:   }
157:   MatDenseRestoreArray(A,&As);
158:   PetscRandomDestroy(&myrand);
159:   MatSetOption(A,MAT_HERMITIAN,PETSC_FALSE);
160:   TestMatInvSqrt(fn,A,viewer,verbose,inplace);

162:   MatDestroy(&A);
163:   FNDestroy(&fn);
164:   SlepcFinalize();
165:   return 0;
166: }

168: /*TEST

170:    testset:
171:       args: -fn_scale 0.9,0.5 -n 10
172:       filter: grep -v "computing matrix functions"
173:       requires: !__float128
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 2 -matcuda
181:          requires: cuda
182:       test:
183:          suffix: 1_magma
184:          args: -fn_method {{1 3}} -matcuda
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 2 -matcuda
192:          requires: cuda
193:       test:
194:          suffix: 2_magma
195:          args: -inplace -fn_method {{1 3}} -matcuda
196:          requires: cuda magma

198: TEST*/