Actual source code: test16.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 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: SlepcInitialize(&argc,&argv,(char*)0,help);
27: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
28: PetscOptionsGetInt(NULL,NULL,"-k",&k,NULL);
29: PetscOptionsGetInt(NULL,NULL,"-l",&l,NULL);
30: PetscOptionsGetInt(NULL,NULL,"-d",&d,NULL);
31: PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
32: PetscPrintf(PETSC_COMM_WORLD,"Test tensor BV of degree %" PetscInt_FMT " with %" PetscInt_FMT " columns of dimension %" PetscInt_FMT "*d.\n",d,k,n);
34: /* Create template vector */
35: VecCreate(PETSC_COMM_WORLD,&t);
36: VecSetSizes(t,PETSC_DECIDE,n);
37: VecSetFromOptions(t);
39: /* Create BV object U */
40: BVCreate(PETSC_COMM_WORLD,&U);
41: PetscObjectSetName((PetscObject)U,"U");
42: BVSetSizesFromVec(U,t,k+d-1);
43: BVSetFromOptions(U);
44: PetscObjectSetName((PetscObject)U,"U");
46: /* Fill first d columns of U */
47: for (j=0;j<d;j++) {
48: BVGetColumn(U,j,&v);
49: VecSet(v,0.0);
50: for (i=0;i<4;i++) {
51: if (i+j<n) VecSetValue(v,i+j,(PetscScalar)(3*i+j-2),INSERT_VALUES);
52: }
53: VecAssemblyBegin(v);
54: VecAssemblyEnd(v);
55: BVRestoreColumn(U,j,&v);
56: }
58: /* Create tensor BV */
59: BVCreateTensor(U,d,&V);
60: PetscObjectSetName((PetscObject)V,"V");
61: BVTensorGetDegree(V,°);
64: /* Set up viewer */
65: PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&view);
66: PetscViewerPushFormat(view,PETSC_VIEWER_ASCII_INFO_DETAIL);
67: BVView(V,view);
68: PetscViewerPopFormat(view);
69: if (verbose) {
70: PetscViewerPushFormat(view,PETSC_VIEWER_ASCII_MATLAB);
71: BVView(V,view);
72: }
74: /* Build first column from previously introduced coefficients */
75: BVTensorBuildFirstColumn(V,d);
76: if (verbose) {
77: PetscPrintf(PETSC_COMM_WORLD,"After building the first column - - - - -\n");
78: BVView(V,view);
79: }
81: /* Test orthogonalization */
82: BVTensorGetFactors(V,&UU,&S);
83: BVGetActiveColumns(UU,NULL,&j);
84: BVGetSizes(UU,NULL,NULL,&id);
86: lds = id*d;
87: for (jj=1;jj<k;jj++) {
88: /* set new orthogonal column in U */
89: BVGetColumn(UU,j,&v);
90: VecSet(v,0.0);
91: for (i=0;i<4;i++) {
92: if (i+j<n) VecSetValue(v,i+j,(PetscScalar)(3*i+j-2),INSERT_VALUES);
93: }
94: VecAssemblyBegin(v);
95: VecAssemblyEnd(v);
96: BVRestoreColumn(UU,j,&v);
97: BVOrthonormalizeColumn(UU,j,PETSC_TRUE,NULL,NULL);
98: j++;
99: BVSetActiveColumns(UU,0,j);
100: /* set new column of S */
101: MatDenseGetArray(S,&pS);
102: for (ii=0;ii<d;ii++) {
103: for (i=0;i<ii+jj+1;i++) {
104: pS[i+ii*id+jj*lds] = (PetscScalar)(2*ii+i+0.5*jj);
105: }
106: }
107: MatDenseRestoreArray(S,&pS);
108: BVOrthonormalizeColumn(V,jj,PETSC_TRUE,NULL,NULL);
109: }
110: BVTensorRestoreFactors(V,&UU,&S);
111: if (verbose) {
112: PetscPrintf(PETSC_COMM_WORLD,"After orthogonalization - - - - -\n");
113: BVView(V,view);
114: }
116: /* Check orthogonality */
117: MatCreateSeqDense(PETSC_COMM_SELF,k,k,NULL,&M);
118: BVDot(V,V,M);
119: MatShift(M,-1.0);
120: MatNorm(M,NORM_1,&norm);
121: if (norm<100*PETSC_MACHINE_EPSILON) PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality < 100*eps\n");
122: else PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality: %g\n",(double)norm);
124: /* Test BVTensorCompress */
125: BVSetActiveColumns(V,0,l);
126: BVTensorCompress(V,0);
127: if (verbose) {
128: PetscPrintf(PETSC_COMM_WORLD,"After BVTensorCompress - - - - -\n");
129: BVView(V,view);
130: }
132: /* Create Mat */
133: MatCreateSeqDense(PETSC_COMM_SELF,k,l,NULL,&Q);
134: PetscObjectSetName((PetscObject)Q,"Q");
135: MatDenseGetArray(Q,&q);
136: for (i=0;i<k;i++)
137: for (j=0;j<l;j++)
138: q[i+j*k] = (i<j)? 2.0: -0.5;
139: MatDenseRestoreArray(Q,&q);
140: if (verbose) MatView(Q,NULL);
142: /* Test BVMultInPlace */
143: BVMultInPlace(V,Q,1,l);
144: if (verbose) {
145: PetscPrintf(PETSC_COMM_WORLD,"After BVMultInPlace - - - - -\n");
146: BVView(V,view);
147: }
149: /* Test BVNorm */
150: BVNorm(V,NORM_1,&norm);
151: PetscPrintf(PETSC_COMM_WORLD,"Norm: %g\n",(double)norm);
153: BVDestroy(&U);
154: BVDestroy(&V);
155: MatDestroy(&Q);
156: MatDestroy(&M);
157: VecDestroy(&t);
158: SlepcFinalize();
159: return 0;
160: }
162: /*TEST
164: test:
165: suffix: 1
166: nsize: 2
167: args: -bv_type {{vecs contiguous svec mat}shared output}
168: output_file: output/test16_1.out
169: filter: grep -v "doing matmult"
171: TEST*/