Actual source code: test4.c

slepc-3.18.2 2023-01-26
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 BV operations, changing the number of active columns.\n\n";

 13: #include <slepcbv.h>

 15: int main(int argc,char **argv)
 16: {
 17:   Vec            t,v;
 18:   Mat            Q,M;
 19:   BV             X,Y;
 20:   PetscInt       i,j,n=10,kx=6,lx=3,ky=5,ly=2;
 21:   PetscScalar    *q,*z;
 22:   PetscReal      nrm;
 23:   PetscViewer    view;
 24:   PetscBool      verbose,trans;

 27:   SlepcInitialize(&argc,&argv,(char*)0,help);
 28:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 29:   PetscOptionsGetInt(NULL,NULL,"-kx",&kx,NULL);
 30:   PetscOptionsGetInt(NULL,NULL,"-lx",&lx,NULL);
 31:   PetscOptionsGetInt(NULL,NULL,"-ky",&ky,NULL);
 32:   PetscOptionsGetInt(NULL,NULL,"-ly",&ly,NULL);
 33:   PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
 34:   PetscPrintf(PETSC_COMM_WORLD,"First BV with %" PetscInt_FMT " active columns (%" PetscInt_FMT " leading columns) of dimension %" PetscInt_FMT ".\n",kx,lx,n);
 35:   PetscPrintf(PETSC_COMM_WORLD,"Second BV with %" PetscInt_FMT " active columns (%" PetscInt_FMT " leading columns) of dimension %" PetscInt_FMT ".\n",ky,ly,n);

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

 42:   /* Create BV object X */
 43:   BVCreate(PETSC_COMM_WORLD,&X);
 44:   PetscObjectSetName((PetscObject)X,"X");
 45:   BVSetSizesFromVec(X,t,kx+2);  /* two extra columns to test active columns */
 46:   BVSetFromOptions(X);
 47:   BVSetActiveColumns(X,lx,kx);

 49:   /* Set up viewer */
 50:   PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&view);
 51:   if (verbose) PetscViewerPushFormat(view,PETSC_VIEWER_ASCII_MATLAB);

 53:   /* Fill X entries */
 54:   for (j=0;j<kx+2;j++) {
 55:     BVGetColumn(X,j,&v);
 56:     VecSet(v,0.0);
 57:     for (i=0;i<4;i++) {
 58:       if (i+j<n) VecSetValue(v,i+j,(PetscScalar)(3*i+j-2),INSERT_VALUES);
 59:     }
 60:     VecAssemblyBegin(v);
 61:     VecAssemblyEnd(v);
 62:     BVRestoreColumn(X,j,&v);
 63:   }
 64:   if (verbose) BVView(X,view);

 66:   /* Create BV object Y */
 67:   BVCreate(PETSC_COMM_WORLD,&Y);
 68:   PetscObjectSetName((PetscObject)Y,"Y");
 69:   BVSetSizesFromVec(Y,t,ky+1);
 70:   BVSetFromOptions(Y);
 71:   BVSetActiveColumns(Y,ly,ky);

 73:   /* Fill Y entries */
 74:   for (j=0;j<ky+1;j++) {
 75:     BVGetColumn(Y,j,&v);
 76:     VecSet(v,(PetscScalar)(j+1)/4.0);
 77:     BVRestoreColumn(Y,j,&v);
 78:   }
 79:   if (verbose) BVView(Y,view);

 81:   /* Create Mat */
 82:   MatCreateSeqDense(PETSC_COMM_SELF,kx,ky,NULL,&Q);
 83:   PetscObjectSetName((PetscObject)Q,"Q");
 84:   MatDenseGetArray(Q,&q);
 85:   for (i=0;i<kx;i++)
 86:     for (j=0;j<ky;j++)
 87:       q[i+j*kx] = (i<j)? 2.0: -0.5;
 88:   MatDenseRestoreArray(Q,&q);
 89:   if (verbose) MatView(Q,NULL);

 91:   /* Test BVResize */
 92:   BVResize(X,kx+4,PETSC_TRUE);

 94:   /* Test BVMult */
 95:   BVMult(Y,2.0,0.5,X,Q);
 96:   if (verbose) {
 97:     PetscPrintf(PETSC_COMM_WORLD,"After BVMult - - - - - - - - -\n");
 98:     BVView(Y,view);
 99:   }

101:   /* Test BVMultVec */
102:   BVGetColumn(Y,0,&v);
103:   PetscMalloc1(kx-lx,&z);
104:   z[0] = 2.0;
105:   for (i=1;i<kx-lx;i++) z[i] = -0.5*z[i-1];
106:   BVMultVec(X,-1.0,1.0,v,z);
107:   PetscFree(z);
108:   BVRestoreColumn(Y,0,&v);
109:   if (verbose) {
110:     PetscPrintf(PETSC_COMM_WORLD,"After BVMultVec - - - - - - -\n");
111:     BVView(Y,view);
112:   }

114:   /* Test BVDot */
115:   MatCreateSeqDense(PETSC_COMM_SELF,ky,kx,NULL,&M);
116:   PetscObjectSetName((PetscObject)M,"M");
117:   BVDot(X,Y,M);
118:   if (verbose) {
119:     PetscPrintf(PETSC_COMM_WORLD,"After BVDot - - - - - - - - -\n");
120:     MatView(M,NULL);
121:   }

123:   /* Test BVDotVec */
124:   BVGetColumn(Y,0,&v);
125:   PetscMalloc1(kx-lx,&z);
126:   BVDotVec(X,v,z);
127:   BVRestoreColumn(Y,0,&v);
128:   if (verbose) {
129:     PetscPrintf(PETSC_COMM_WORLD,"After BVDotVec - - - - - - -\n");
130:     VecCreateSeqWithArray(PETSC_COMM_SELF,1,kx-lx,z,&v);
131:     PetscObjectSetName((PetscObject)v,"z");
132:     VecView(v,view);
133:     VecDestroy(&v);
134:   }
135:   PetscFree(z);

137:   /* Test BVMultInPlace and BVScale */
138:   PetscOptionsHasName(NULL,NULL,"-trans",&trans);
139:   if (trans) {
140:     Mat Qt;
141:     MatTranspose(Q,MAT_INITIAL_MATRIX,&Qt);
142:     BVMultInPlaceHermitianTranspose(X,Qt,lx+1,ky);
143:     MatDestroy(&Qt);
144:   } else BVMultInPlace(X,Q,lx+1,ky);
145:   BVScale(X,2.0);
146:   if (verbose) {
147:     PetscPrintf(PETSC_COMM_WORLD,"After BVMultInPlace - - - - -\n");
148:     BVView(X,view);
149:   }

151:   /* Test BVNorm */
152:   BVNormColumn(X,lx,NORM_2,&nrm);
153:   PetscPrintf(PETSC_COMM_WORLD,"2-Norm of X[%" PetscInt_FMT "] = %g\n",lx,(double)nrm);
154:   BVNorm(X,NORM_FROBENIUS,&nrm);
155:   PetscPrintf(PETSC_COMM_WORLD,"Frobenius Norm of X = %g\n",(double)nrm);

157:   BVDestroy(&X);
158:   BVDestroy(&Y);
159:   MatDestroy(&Q);
160:   MatDestroy(&M);
161:   VecDestroy(&t);
162:   SlepcFinalize();
163:   return 0;
164: }

166: /*TEST

168:    testset:
169:       output_file: output/test4_1.out
170:       args: -n 18 -kx 12 -ky 8
171:       test:
172:          suffix: 1
173:          args: -bv_type {{vecs contiguous svec mat}shared output}
174:       test:
175:          suffix: 1_vecs_vmip
176:          args: -bv_type vecs -bv_vecs_vmip 1
177:       test:
178:          suffix: 1_cuda
179:          args: -bv_type svec -vec_type cuda
180:          requires: cuda
181:       test:
182:          suffix: 2
183:          args: -bv_type {{vecs contiguous svec mat}shared output} -trans

185: TEST*/