Actual source code: test4.c

slepc-3.22.1 2024-10-28
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;

 26:   PetscFunctionBeginUser;
 27:   PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
 28:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
 29:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-kx",&kx,NULL));
 30:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-lx",&lx,NULL));
 31:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-ky",&ky,NULL));
 32:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-ly",&ly,NULL));
 33:   PetscCall(PetscOptionsHasName(NULL,NULL,"-verbose",&verbose));
 34:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"First BV with %" PetscInt_FMT " active columns (%" PetscInt_FMT " leading columns) of dimension %" PetscInt_FMT ".\n",kx,lx,n));
 35:   PetscCall(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:   PetscCall(VecCreate(PETSC_COMM_WORLD,&t));
 39:   PetscCall(VecSetSizes(t,PETSC_DECIDE,n));
 40:   PetscCall(VecSetFromOptions(t));

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

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

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

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

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

 81:   /* Create Mat */
 82:   PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,kx,ky,NULL,&Q));
 83:   PetscCall(PetscObjectSetName((PetscObject)Q,"Q"));
 84:   PetscCall(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:   PetscCall(MatDenseRestoreArray(Q,&q));
 89:   if (verbose) PetscCall(MatView(Q,NULL));

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

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

101:   /* Test BVMultVec */
102:   PetscCall(BVGetColumn(Y,0,&v));
103:   PetscCall(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:   PetscCall(BVMultVec(X,-1.0,1.0,v,z));
107:   PetscCall(PetscFree(z));
108:   PetscCall(BVRestoreColumn(Y,0,&v));
109:   if (verbose) {
110:     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"After BVMultVec - - - - - - -\n"));
111:     PetscCall(BVView(Y,view));
112:   }

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

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

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

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

157:   PetscCall(BVDestroy(&X));
158:   PetscCall(BVDestroy(&Y));
159:   PetscCall(MatDestroy(&Q));
160:   PetscCall(MatDestroy(&M));
161:   PetscCall(VecDestroy(&t));
162:   PetscCall(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 mat}} -vec_type cuda
180:          requires: cuda
181:       test:
182:          suffix: 1_hip
183:          args: -bv_type {{svec mat}} -vec_type hip
184:          requires: hip
185:       test:
186:          suffix: 2
187:          args: -bv_type {{vecs contiguous svec mat}shared output} -trans

189: TEST*/