Actual source code: test16.c

slepc-3.21.2 2024-09-25
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 tensor BV.\n\n";

 13: #include <slepcbv.h>

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

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

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

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

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

 59:   /* Create tensor BV */
 60:   PetscCall(BVCreateTensor(U,d,&V));
 61:   PetscCall(PetscObjectSetName((PetscObject)V,"V"));
 62:   PetscCall(BVTensorGetDegree(V,&deg));
 63:   PetscCheck(deg==d,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Wrong degree");

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

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

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

117:   /* Check orthogonality */
118:   PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,k,k,NULL,&M));
119:   PetscCall(BVDot(V,V,M));
120:   PetscCall(MatShift(M,-1.0));
121:   PetscCall(MatNorm(M,NORM_1,&norm));
122:   if (norm<100*PETSC_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality < 100*eps\n"));
123:   else PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality: %g\n",(double)norm));

125:   /* Test BVTensorCompress */
126:   PetscCall(BVSetActiveColumns(V,0,l));
127:   PetscCall(BVTensorCompress(V,0));
128:   if (verbose) {
129:     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"After BVTensorCompress - - - - -\n"));
130:     PetscCall(BVView(V,view));
131:   }

133:   /* Create Mat */
134:   PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,k,l,NULL,&Q));
135:   PetscCall(PetscObjectSetName((PetscObject)Q,"Q"));
136:   PetscCall(MatDenseGetArray(Q,&q));
137:   for (i=0;i<k;i++)
138:     for (j=0;j<l;j++)
139:       q[i+j*k] = (i<j)? 2.0: -0.5;
140:   PetscCall(MatDenseRestoreArray(Q,&q));
141:   if (verbose) PetscCall(MatView(Q,NULL));

143:   /* Test BVMultInPlace */
144:   PetscCall(BVMultInPlace(V,Q,1,l));
145:   if (verbose) {
146:     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"After BVMultInPlace - - - - -\n"));
147:     PetscCall(BVView(V,view));
148:   }

150:   /* Test BVNorm */
151:   PetscCall(BVNorm(V,NORM_1,&norm));
152:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Norm: %g\n",(double)norm));

154:   PetscCall(BVDestroy(&U));
155:   PetscCall(BVDestroy(&V));
156:   PetscCall(MatDestroy(&Q));
157:   PetscCall(MatDestroy(&M));
158:   PetscCall(VecDestroy(&t));
159:   PetscCall(SlepcFinalize());
160:   return 0;
161: }

163: /*TEST

165:    testset:
166:       nsize: 2
167:       output_file: output/test16_1.out
168:       filter: grep -v "doing matmult"
169:       test:
170:          suffix: 1
171:          args: -bv_type {{vecs contiguous svec mat}}
172:       test:
173:          suffix: 1_cuda
174:          args: -bv_type {{vecs svec mat}} -vec_type cuda
175:          requires: cuda

177: TEST*/