Actual source code: test13.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 logarithm.\n\n";

 13: #include <slepcfn.h>

 15: /*
 16:    Compute matrix logarithm B = logm(A)
 17:  */
 18: PetscErrorCode TestMatLog(FN fn,Mat A,PetscViewer viewer,PetscBool verbose,PetscBool inplace)
 19: {
 21:   PetscBool      set,flg;
 22:   PetscScalar    tau,eta;
 23:   PetscInt       n;
 24:   Mat            F,R;
 25:   Vec            v,f0;
 26:   FN             fnexp;
 27:   PetscReal      nrm;

 30:   MatGetSize(A,&n,NULL);
 31:   MatCreateSeqDense(PETSC_COMM_SELF,n,n,NULL,&F);
 32:   PetscObjectSetName((PetscObject)F,"F");
 33:   MatCreateSeqDense(PETSC_COMM_SELF,n,n,NULL,&R);
 34:   PetscObjectSetName((PetscObject)R,"R");
 35:   FNGetScale(fn,&tau,&eta);
 36:   /* compute matrix logarithm */
 37:   if (inplace) {
 38:     MatCopy(A,F,SAME_NONZERO_PATTERN);
 39:     MatIsHermitianKnown(A,&set,&flg);
 40:     if (set && flg) { MatSetOption(F,MAT_HERMITIAN,PETSC_TRUE); }
 41:     FNEvaluateFunctionMat(fn,F,NULL);
 42:   } else {
 43:     FNEvaluateFunctionMat(fn,A,F);
 44:   }
 45:   if (verbose) {
 46:     PetscPrintf(PETSC_COMM_WORLD,"Matrix A - - - - - - - -\n");
 47:     MatView(A,viewer);
 48:     PetscPrintf(PETSC_COMM_WORLD,"Computed logm(A) - - - - - - -\n");
 49:     MatView(F,viewer);
 50:   }
 51:   /* check error ||expm(F)-A||_F */
 52:   FNCreate(PETSC_COMM_WORLD,&fnexp);
 53:   FNSetType(fnexp,FNEXP);
 54:   MatCopy(F,R,SAME_NONZERO_PATTERN);
 55:   if (eta!=1.0) {
 56:     MatScale(R,1.0/eta);
 57:   }
 58:   FNEvaluateFunctionMat(fnexp,R,NULL);
 59:   FNDestroy(&fnexp);
 60:   MatAXPY(R,-tau,A,SAME_NONZERO_PATTERN);
 61:   MatNorm(R,NORM_FROBENIUS,&nrm);
 62:   if (nrm<100*PETSC_MACHINE_EPSILON) {
 63:     PetscPrintf(PETSC_COMM_WORLD,"||expm(F)-A||_F < 100*eps\n");
 64:   } else {
 65:     PetscPrintf(PETSC_COMM_WORLD,"||expm(F)-A||_F = %g\n",(double)nrm);
 66:   }
 67:   /* check FNEvaluateFunctionMatVec() */
 68:   MatCreateVecs(A,&v,&f0);
 69:   MatGetColumnVector(F,f0,0);
 70:   FNEvaluateFunctionMatVec(fn,A,v);
 71:   VecAXPY(v,-1.0,f0);
 72:   VecNorm(v,NORM_2,&nrm);
 73:   if (nrm>100*PETSC_MACHINE_EPSILON) {
 74:     PetscPrintf(PETSC_COMM_WORLD,"Warning: the norm of f(A)*e_1-v is %g\n",(double)nrm);
 75:   }
 76:   MatDestroy(&F);
 77:   MatDestroy(&R);
 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    *As;
 90:   PetscViewer    viewer;
 91:   PetscBool      verbose,inplace,random;

 93:   SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
 94:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 95:   PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
 96:   PetscOptionsHasName(NULL,NULL,"-inplace",&inplace);
 97:   PetscOptionsHasName(NULL,NULL,"-random",&random);
 98:   PetscPrintf(PETSC_COMM_WORLD,"Matrix logarithm, n=%D.\n",n);

100:   /* Create logarithm function object */
101:   FNCreate(PETSC_COMM_WORLD,&fn);
102:   FNSetType(fn,FNLOG);
103:   FNSetFromOptions(fn);

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

112:   /* Create matrices */
113:   MatCreateSeqDense(PETSC_COMM_SELF,n,n,NULL,&A);
114:   PetscObjectSetName((PetscObject)A,"A");

116:   if (random) {
117:     MatSetRandom(A,NULL);
118:   } else {
119:     /* Fill A with a non-symmetric Toeplitz matrix */
120:     MatDenseGetArray(A,&As);
121:     for (i=0;i<n;i++) As[i+i*n]=2.0;
122:     for (j=1;j<3;j++) {
123:       for (i=0;i<n-j;i++) { As[i+(i+j)*n]=1.0; As[(i+j)+i*n]=-1.0; }
124:     }
125:     As[(n-1)*n] = -5.0;
126:     MatDenseRestoreArray(A,&As);
127:   }
128:   TestMatLog(fn,A,viewer,verbose,inplace);

130:   MatDestroy(&A);
131:   FNDestroy(&fn);
132:   SlepcFinalize();
133:   return ierr;
134: }

136: /*TEST

138:    testset:
139:       filter: grep -v "computing matrix functions"
140:       output_file: output/test13_1.out
141:       test:
142:          suffix: 1
143:          args: -fn_scale .04,2 -n 75
144:          requires: c99_complex !__float128
145:       test:
146:          suffix: 1_random
147:          args: -fn_scale .04,2 -n 75 -random
148:          requires: complex

150: TEST*/