Actual source code: test16.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 tensor BV.\n\n";

 13: #include <slepcbv.h>

 15: int main(int argc,char **argv)
 16: {
 17:   PetscErrorCode    ierr;
 18:   Vec               t,v;
 19:   Mat               S,M,Q;
 20:   BV                U,V,UU;
 21:   PetscInt          i,ii,j,jj,n=10,k=6,l=3,d=3,deg,id,lds;
 22:   PetscScalar       *pS,*q;
 23:   PetscReal         norm;
 24:   PetscViewer       view;
 25:   PetscBool         verbose;

 27:   SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
 28:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 29:   PetscOptionsGetInt(NULL,NULL,"-k",&k,NULL);
 30:   PetscOptionsGetInt(NULL,NULL,"-l",&l,NULL);
 31:   PetscOptionsGetInt(NULL,NULL,"-d",&d,NULL);
 32:   PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
 33:   PetscPrintf(PETSC_COMM_WORLD,"Test tensor BV of degree %D with %D columns of dimension %D*d.\n",d,k,n);

 35:   /* Create template vector */
 36:   VecCreate(PETSC_COMM_WORLD,&t);
 37:   VecSetSizes(t,PETSC_DECIDE,n);
 38:   VecSetFromOptions(t);

 40:   /* Create BV object U */
 41:   BVCreate(PETSC_COMM_WORLD,&U);
 42:   PetscObjectSetName((PetscObject)U,"U");
 43:   BVSetSizesFromVec(U,t,k+d-1);
 44:   BVSetFromOptions(U);
 45:   PetscObjectSetName((PetscObject)U,"U");

 47:   /* Fill first d columns of U */
 48:   for (j=0;j<d;j++) {
 49:     BVGetColumn(U,j,&v);
 50:     VecSet(v,0.0);
 51:     for (i=0;i<4;i++) {
 52:       if (i+j<n) {
 53:         VecSetValue(v,i+j,(PetscScalar)(3*i+j-2),INSERT_VALUES);
 54:       }
 55:     }
 56:     VecAssemblyBegin(v);
 57:     VecAssemblyEnd(v);
 58:     BVRestoreColumn(U,j,&v);
 59:   }

 61:   /* Create tensor BV */
 62:   BVCreateTensor(U,d,&V);
 63:   PetscObjectSetName((PetscObject)V,"V");
 64:   BVTensorGetDegree(V,&deg);
 65:   if (deg!=d) SETERRQ(PETSC_COMM_WORLD,1,"Wrong degree");

 67:   /* Set up viewer */
 68:   PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&view);
 69:   PetscViewerPushFormat(view,PETSC_VIEWER_ASCII_INFO_DETAIL);
 70:   BVView(V,view);
 71:   PetscViewerPopFormat(view);
 72:   if (verbose) {
 73:     PetscViewerPushFormat(view,PETSC_VIEWER_ASCII_MATLAB);
 74:     BVView(V,view);
 75:   }

 77:   /* Build first column from previously introduced coefficients */
 78:   BVTensorBuildFirstColumn(V,d);
 79:   if (verbose) {
 80:     PetscPrintf(PETSC_COMM_WORLD,"After building the first column - - - - -\n");
 81:     BVView(V,view);
 82:   }

 84:   /* Test orthogonalization */
 85:   BVTensorGetFactors(V,&UU,&S);
 86:   BVGetActiveColumns(UU,NULL,&j);
 87:   BVGetSizes(UU,NULL,NULL,&id);
 88:   if (id!=k+d-1) SETERRQ(PETSC_COMM_WORLD,1,"Wrong dimensions");
 89:   lds = id*d;
 90:   for (jj=1;jj<k;jj++) {
 91:     /* set new orthogonal column in U */
 92:     BVGetColumn(UU,j,&v);
 93:     VecSet(v,0.0);
 94:     for (i=0;i<4;i++) {
 95:       if (i+j<n) {
 96:         VecSetValue(v,i+j,(PetscScalar)(3*i+j-2),INSERT_VALUES);
 97:       }
 98:     }
 99:     VecAssemblyBegin(v);
100:     VecAssemblyEnd(v);
101:     BVRestoreColumn(UU,j,&v);
102:     BVOrthonormalizeColumn(UU,j,PETSC_TRUE,NULL,NULL);
103:     j++;
104:     BVSetActiveColumns(UU,0,j);
105:     /* set new column of S */
106:     MatDenseGetArray(S,&pS);
107:     for (ii=0;ii<d;ii++) {
108:       for (i=0;i<ii+jj+1;i++) {
109:         pS[i+ii*id+jj*lds] = (PetscScalar)(2*ii+i+0.5*jj);
110:       }
111:     }
112:     MatDenseRestoreArray(S,&pS);
113:     BVOrthonormalizeColumn(V,jj,PETSC_TRUE,NULL,NULL);
114:   }
115:   BVTensorRestoreFactors(V,&UU,&S);
116:   if (verbose) {
117:     PetscPrintf(PETSC_COMM_WORLD,"After orthogonalization - - - - -\n");
118:     BVView(V,view);
119:   }

121:   /* Check orthogonality */
122:   MatCreateSeqDense(PETSC_COMM_SELF,k,k,NULL,&M);
123:   BVDot(V,V,M);
124:   MatShift(M,-1.0);
125:   MatNorm(M,NORM_1,&norm);
126:   if (norm<100*PETSC_MACHINE_EPSILON) {
127:     PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality < 100*eps\n");
128:   } else {
129:     PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality: %g\n",(double)norm);
130:   }

132:   /* Test BVTensorCompress */
133:   BVSetActiveColumns(V,0,l);
134:   BVTensorCompress(V,0);
135:   if (verbose) {
136:     PetscPrintf(PETSC_COMM_WORLD,"After BVTensorCompress - - - - -\n");
137:     BVView(V,view);
138:   }

140:   /* Create Mat */
141:   MatCreateSeqDense(PETSC_COMM_SELF,k,l,NULL,&Q);
142:   PetscObjectSetName((PetscObject)Q,"Q");
143:   MatDenseGetArray(Q,&q);
144:   for (i=0;i<k;i++)
145:     for (j=0;j<l;j++)
146:       q[i+j*k] = (i<j)? 2.0: -0.5;
147:   MatDenseRestoreArray(Q,&q);
148:   if (verbose) {
149:     MatView(Q,NULL);
150:   }

152:   /* Test BVMultInPlace */
153:   BVMultInPlace(V,Q,1,l);
154:   if (verbose) {
155:     PetscPrintf(PETSC_COMM_WORLD,"After BVMultInPlace - - - - -\n");
156:     BVView(V,view);
157:   }

159:   BVDestroy(&U);
160:   BVDestroy(&V);
161:   MatDestroy(&Q);
162:   MatDestroy(&M);
163:   VecDestroy(&t);
164:   SlepcFinalize();
165:   return ierr;
166: }

168: /*TEST

170:    test:
171:       suffix: 1
172:       nsize: 2
173:       args: -bv_type {{vecs contiguous svec mat}shared output}
174:       output_file: output/test16_1.out

176: TEST*/