Actual source code: test4.c
slepc-3.17.1 2022-04-11
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: SlepcInitialize(&argc,&argv,(char*)0,help);
27: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
28: PetscOptionsGetInt(NULL,NULL,"-kx",&kx,NULL);
29: PetscOptionsGetInt(NULL,NULL,"-lx",&lx,NULL);
30: PetscOptionsGetInt(NULL,NULL,"-ky",&ky,NULL);
31: PetscOptionsGetInt(NULL,NULL,"-ly",&ly,NULL);
32: PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
33: PetscPrintf(PETSC_COMM_WORLD,"First BV with %" PetscInt_FMT " active columns (%" PetscInt_FMT " leading columns) of dimension %" PetscInt_FMT ".\n",kx,lx,n);
34: PetscPrintf(PETSC_COMM_WORLD,"Second BV with %" PetscInt_FMT " active columns (%" PetscInt_FMT " leading columns) of dimension %" PetscInt_FMT ".\n",ky,ly,n);
36: /* Create template vector */
37: VecCreate(PETSC_COMM_WORLD,&t);
38: VecSetSizes(t,PETSC_DECIDE,n);
39: VecSetFromOptions(t);
41: /* Create BV object X */
42: BVCreate(PETSC_COMM_WORLD,&X);
43: PetscObjectSetName((PetscObject)X,"X");
44: BVSetSizesFromVec(X,t,kx+2); /* two extra columns to test active columns */
45: BVSetFromOptions(X);
46: BVSetActiveColumns(X,lx,kx);
48: /* Set up viewer */
49: PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&view);
50: if (verbose) PetscViewerPushFormat(view,PETSC_VIEWER_ASCII_MATLAB);
52: /* Fill X entries */
53: for (j=0;j<kx+2;j++) {
54: BVGetColumn(X,j,&v);
55: VecSet(v,0.0);
56: for (i=0;i<4;i++) {
57: if (i+j<n) VecSetValue(v,i+j,(PetscScalar)(3*i+j-2),INSERT_VALUES);
58: }
59: VecAssemblyBegin(v);
60: VecAssemblyEnd(v);
61: BVRestoreColumn(X,j,&v);
62: }
63: if (verbose) BVView(X,view);
65: /* Create BV object Y */
66: BVCreate(PETSC_COMM_WORLD,&Y);
67: PetscObjectSetName((PetscObject)Y,"Y");
68: BVSetSizesFromVec(Y,t,ky+1);
69: BVSetFromOptions(Y);
70: BVSetActiveColumns(Y,ly,ky);
72: /* Fill Y entries */
73: for (j=0;j<ky+1;j++) {
74: BVGetColumn(Y,j,&v);
75: VecSet(v,(PetscScalar)(j+1)/4.0);
76: BVRestoreColumn(Y,j,&v);
77: }
78: if (verbose) BVView(Y,view);
80: /* Create Mat */
81: MatCreateSeqDense(PETSC_COMM_SELF,kx,ky,NULL,&Q);
82: PetscObjectSetName((PetscObject)Q,"Q");
83: MatDenseGetArray(Q,&q);
84: for (i=0;i<kx;i++)
85: for (j=0;j<ky;j++)
86: q[i+j*kx] = (i<j)? 2.0: -0.5;
87: MatDenseRestoreArray(Q,&q);
88: if (verbose) MatView(Q,NULL);
90: /* Test BVResize */
91: BVResize(X,kx+4,PETSC_TRUE);
93: /* Test BVMult */
94: BVMult(Y,2.0,0.5,X,Q);
95: if (verbose) {
96: PetscPrintf(PETSC_COMM_WORLD,"After BVMult - - - - - - - - -\n");
97: BVView(Y,view);
98: }
100: /* Test BVMultVec */
101: BVGetColumn(Y,0,&v);
102: PetscMalloc1(kx-lx,&z);
103: z[0] = 2.0;
104: for (i=1;i<kx-lx;i++) z[i] = -0.5*z[i-1];
105: BVMultVec(X,-1.0,1.0,v,z);
106: PetscFree(z);
107: BVRestoreColumn(Y,0,&v);
108: if (verbose) {
109: PetscPrintf(PETSC_COMM_WORLD,"After BVMultVec - - - - - - -\n");
110: BVView(Y,view);
111: }
113: /* Test BVDot */
114: MatCreateSeqDense(PETSC_COMM_SELF,ky,kx,NULL,&M);
115: PetscObjectSetName((PetscObject)M,"M");
116: BVDot(X,Y,M);
117: if (verbose) {
118: PetscPrintf(PETSC_COMM_WORLD,"After BVDot - - - - - - - - -\n");
119: MatView(M,NULL);
120: }
122: /* Test BVDotVec */
123: BVGetColumn(Y,0,&v);
124: PetscMalloc1(kx-lx,&z);
125: BVDotVec(X,v,z);
126: BVRestoreColumn(Y,0,&v);
127: if (verbose) {
128: PetscPrintf(PETSC_COMM_WORLD,"After BVDotVec - - - - - - -\n");
129: VecCreateSeqWithArray(PETSC_COMM_SELF,1,kx-lx,z,&v);
130: PetscObjectSetName((PetscObject)v,"z");
131: VecView(v,view);
132: VecDestroy(&v);
133: }
134: PetscFree(z);
136: /* Test BVMultInPlace and BVScale */
137: PetscOptionsHasName(NULL,NULL,"-trans",&trans);
138: if (trans) {
139: Mat Qt;
140: MatTranspose(Q,MAT_INITIAL_MATRIX,&Qt);
141: BVMultInPlaceHermitianTranspose(X,Qt,lx+1,ky);
142: MatDestroy(&Qt);
143: } else BVMultInPlace(X,Q,lx+1,ky);
144: BVScale(X,2.0);
145: if (verbose) {
146: PetscPrintf(PETSC_COMM_WORLD,"After BVMultInPlace - - - - -\n");
147: BVView(X,view);
148: }
150: /* Test BVNorm */
151: BVNormColumn(X,lx,NORM_2,&nrm);
152: PetscPrintf(PETSC_COMM_WORLD,"2-Norm of X[%" PetscInt_FMT "] = %g\n",lx,(double)nrm);
153: BVNorm(X,NORM_FROBENIUS,&nrm);
154: PetscPrintf(PETSC_COMM_WORLD,"Frobenius Norm of X = %g\n",(double)nrm);
156: BVDestroy(&X);
157: BVDestroy(&Y);
158: MatDestroy(&Q);
159: MatDestroy(&M);
160: VecDestroy(&t);
161: SlepcFinalize();
162: return 0;
163: }
165: /*TEST
167: testset:
168: output_file: output/test4_1.out
169: args: -n 18 -kx 12 -ky 8
170: test:
171: suffix: 1
172: args: -bv_type {{vecs contiguous svec mat}shared output}
173: test:
174: suffix: 1_vecs_vmip
175: args: -bv_type vecs -bv_vecs_vmip 1
176: test:
177: suffix: 1_cuda
178: args: -bv_type svec -vec_type cuda
179: requires: cuda
180: test:
181: suffix: 2
182: args: -bv_type {{vecs contiguous svec mat}shared output} -trans
184: TEST*/