Actual source code: test8.c

slepc-3.7.2 2016-07-19
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2016, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.

  8:    SLEPc is free software: you can redistribute it and/or modify it under  the
  9:    terms of version 3 of the GNU Lesser General Public License as published by
 10:    the Free Software Foundation.

 12:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
 13:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
 14:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
 15:    more details.

 17:    You  should have received a copy of the GNU Lesser General  Public  License
 18:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 19:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 20: */

 22: static char help[] = "Test matrix inverse square root.\n\n";

 24: #include <slepcfn.h>

 28: /*
 29:    Compute matrix inverse square root B = inv(sqrtm(A))
 30:    Check result as norm(B*B*A-I)
 31:  */
 32: PetscErrorCode TestMatInvSqrt(FN fn,Mat A,PetscViewer viewer,PetscBool verbose,PetscBool inplace)
 33: {
 35:   PetscScalar    tau,eta;
 36:   PetscReal      nrm;
 37:   PetscBool      set,flg;
 38:   PetscInt       n;
 39:   Mat            S,R;
 40:   Vec            v,f0;

 43:   MatGetSize(A,&n,NULL);
 44:   MatCreateSeqDense(PETSC_COMM_SELF,n,n,NULL,&S);
 45:   PetscObjectSetName((PetscObject)S,"S");
 46:   MatCreateSeqDense(PETSC_COMM_SELF,n,n,NULL,&R);
 47:   PetscObjectSetName((PetscObject)R,"R");
 48:   FNGetScale(fn,&tau,&eta);
 49:   /* compute inverse square root */
 50:   if (inplace) {
 51:     MatCopy(A,S,SAME_NONZERO_PATTERN);
 52:     MatIsHermitianKnown(A,&set,&flg);
 53:     if (set && flg) { MatSetOption(S,MAT_HERMITIAN,PETSC_TRUE); }
 54:     FNEvaluateFunctionMat(fn,S,NULL);
 55:   } else {
 56:     FNEvaluateFunctionMat(fn,A,S);
 57:   }
 58:   if (verbose) {
 59:     PetscPrintf(PETSC_COMM_WORLD,"Matrix A - - - - - - - -\n");
 60:     MatView(A,viewer);
 61:     PetscPrintf(PETSC_COMM_WORLD,"Computed inv(sqrtm(A)) - - - - - - -\n");
 62:     MatView(S,viewer);
 63:   }
 64:   /* check error ||S*S*A-I||_F */
 65:   MatMatMult(S,S,MAT_REUSE_MATRIX,PETSC_DEFAULT,&R);
 66:   if (eta!=1.0) {
 67:     MatScale(R,1.0/(eta*eta));
 68:   }
 69:   MatCreateVecs(A,&v,&f0);
 70:   MatGetColumnVector(S,f0,0);
 71:   MatCopy(R,S,SAME_NONZERO_PATTERN);
 72:   if (tau!=1.0) {
 73:     MatScale(S,tau);
 74:   }
 75:   MatMatMult(S,A,MAT_REUSE_MATRIX,PETSC_DEFAULT,&R);
 76:   MatShift(R,-1.0);
 77:   MatNorm(R,NORM_FROBENIUS,&nrm);
 78:   if (nrm<100*PETSC_MACHINE_EPSILON) {
 79:     PetscPrintf(PETSC_COMM_WORLD,"||S*S*A-I||_F < 100*eps\n");
 80:   } else {
 81:     PetscPrintf(PETSC_COMM_WORLD,"||S*S*A-I||_F = %g\n",(double)nrm);
 82:   }
 83:   /* check FNEvaluateFunctionMatVec() */
 84:   FNEvaluateFunctionMatVec(fn,A,v);
 85:   VecAXPY(v,-1.0,f0);
 86:   VecNorm(v,NORM_2,&nrm);
 87:   if (nrm>100*PETSC_MACHINE_EPSILON) {
 88:     PetscPrintf(PETSC_COMM_WORLD,"Warning: the norm of f(A)*e_1-v is %g\n",(double)nrm);
 89:   }
 90:   MatDestroy(&S);
 91:   MatDestroy(&R);
 92:   VecDestroy(&v);
 93:   VecDestroy(&f0);
 94:   return(0);
 95: }

 99: int main(int argc,char **argv)
100: {
102:   FN             fn;
103:   Mat            A;
104:   PetscInt       i,j,n=10;
105:   PetscScalar    *As,tau=1.0,eta=1.0;
106:   PetscViewer    viewer;
107:   PetscBool      verbose,inplace;
108:   PetscRandom    myrand;
109:   PetscReal      v;

111:   SlepcInitialize(&argc,&argv,(char*)0,help);
112:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
113:   PetscOptionsGetScalar(NULL,NULL,"-tau",&tau,NULL);
114:   PetscOptionsGetScalar(NULL,NULL,"-eta",&eta,NULL);
115:   PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
116:   PetscOptionsHasName(NULL,NULL,"-inplace",&inplace);
117:   PetscPrintf(PETSC_COMM_WORLD,"Matrix inverse square root, n=%D.\n",n);

119:   /* Create function eta*inv(sqrt(tau*x)) */
120:   FNCreate(PETSC_COMM_WORLD,&fn);
121:   FNSetType(fn,FNINVSQRT);
122:   FNSetScale(fn,tau,eta);

124:   /* Set up viewer */
125:   PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
126:   FNView(fn,viewer);
127:   if (verbose) {
128:     PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);
129:   }

131:   /* Create matrix */
132:   MatCreateSeqDense(PETSC_COMM_SELF,n,n,NULL,&A);
133:   PetscObjectSetName((PetscObject)A,"A");

135:   /* Compute square root of a symmetric matrix A */
136:   MatDenseGetArray(A,&As);
137:   for (i=0;i<n;i++) As[i+i*n]=2.5;
138:   for (j=1;j<3;j++) {
139:     for (i=0;i<n-j;i++) { As[i+(i+j)*n]=1.0; As[(i+j)+i*n]=1.0; }
140:   }
141:   MatDenseRestoreArray(A,&As);
142:   MatSetOption(A,MAT_HERMITIAN,PETSC_TRUE);
143:   TestMatInvSqrt(fn,A,viewer,verbose,inplace);

145:   /* Repeat with upper triangular A */
146:   MatDenseGetArray(A,&As);
147:   for (j=1;j<3;j++) {
148:     for (i=0;i<n-j;i++) As[(i+j)+i*n]=0.0;
149:   }
150:   MatDenseRestoreArray(A,&As);
151:   MatSetOption(A,MAT_HERMITIAN,PETSC_FALSE);
152:   TestMatInvSqrt(fn,A,viewer,verbose,inplace);

154:   /* Repeat with non-symmetic A */
155:   PetscRandomCreate(PETSC_COMM_WORLD,&myrand);
156:   PetscRandomSetFromOptions(myrand);
157:   PetscRandomSetInterval(myrand,0.0,1.0);
158:   MatDenseGetArray(A,&As);
159:   for (j=1;j<3;j++) {
160:     for (i=0;i<n-j;i++) { 
161:       PetscRandomGetValueReal(myrand,&v);
162:       As[(i+j)+i*n]=v;
163:     }
164:   }
165:   MatDenseRestoreArray(A,&As);
166:   PetscRandomDestroy(&myrand);
167:   MatSetOption(A,MAT_HERMITIAN,PETSC_FALSE);
168:   TestMatInvSqrt(fn,A,viewer,verbose,inplace);

170:   MatDestroy(&A);
171:   FNDestroy(&fn);
172:   SlepcFinalize();
173:   return ierr;
174: }