Actual source code: test7.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 square root.\n\n";

 24: #include <slepcfn.h>

 28: /*
 29:    Compute matrix square root B = sqrtm(A)
 30:    Check result as norm(B*B-A)
 31:  */
 32: PetscErrorCode TestMatSqrt(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 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 sqrtm(A) - - - - - - -\n");
 62:     MatView(S,viewer);
 63:   }
 64:   /* check error ||S*S-A||_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:   MatAXPY(R,-tau,A,SAME_NONZERO_PATTERN);
 70:   MatNorm(R,NORM_FROBENIUS,&nrm);
 71:   if (nrm<100*PETSC_MACHINE_EPSILON) {
 72:     PetscPrintf(PETSC_COMM_WORLD,"||S*S-A||_F < 100*eps\n");
 73:   } else {
 74:     PetscPrintf(PETSC_COMM_WORLD,"||S*S-A||_F = %g\n",(double)nrm);
 75:   }
 76:   /* check FNEvaluateFunctionMatVec() */
 77:   MatCreateVecs(A,&v,&f0);
 78:   MatGetColumnVector(S,f0,0);
 79:   FNEvaluateFunctionMatVec(fn,A,v);
 80:   VecAXPY(v,-1.0,f0);
 81:   VecNorm(v,NORM_2,&nrm);
 82:   if (nrm>100*PETSC_MACHINE_EPSILON) {
 83:     PetscPrintf(PETSC_COMM_WORLD,"Warning: the norm of f(A)*e_1-v is %g\n",(double)nrm);
 84:   }
 85:   MatDestroy(&S);
 86:   MatDestroy(&R);
 87:   VecDestroy(&v);
 88:   VecDestroy(&f0);
 89:   return(0);
 90: }

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

106:   SlepcInitialize(&argc,&argv,(char*)0,help);
107:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
108:   PetscOptionsGetScalar(NULL,NULL,"-tau",&tau,NULL);
109:   PetscOptionsGetScalar(NULL,NULL,"-eta",&eta,NULL);
110:   PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
111:   PetscOptionsHasName(NULL,NULL,"-inplace",&inplace);
112:   PetscPrintf(PETSC_COMM_WORLD,"Matrix square root, n=%D.\n",n);

114:   /* Create function eta*sqrt(tau*x) */
115:   FNCreate(PETSC_COMM_WORLD,&fn);
116:   FNSetType(fn,FNSQRT);
117:   FNSetScale(fn,tau,eta);

119:   /* Set up viewer */
120:   PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
121:   FNView(fn,viewer);
122:   if (verbose) {
123:     PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);
124:   }

126:   /* Create matrix */
127:   MatCreateSeqDense(PETSC_COMM_SELF,n,n,NULL,&A);
128:   PetscObjectSetName((PetscObject)A,"A");

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

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

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

165:   MatDestroy(&A);
166:   FNDestroy(&fn);
167:   SlepcFinalize();
168:   return ierr;
169: }