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