Actual source code: test19.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 BVGetSplitRows().\n\n";

 13: #include <slepcbv.h>

 15: int main(int argc,char **argv)
 16: {
 17:   Mat            M,T,D,block[4];
 18:   Vec            t,v,v1,v2,w,*C;
 19:   BV             X,U,L;
 20:   IS             is[2];
 21:   PetscInt       i,j,n=10,k=5,l=3,nc=0,rstart,rend;
 22:   PetscViewer    view;
 23:   PetscBool      verbose;

 25:   PetscFunctionBeginUser;
 26:   PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
 27:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
 28:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-k",&k,NULL));
 29:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-l",&l,NULL));
 30:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-nc",&nc,NULL));
 31:   PetscCall(PetscOptionsHasName(NULL,NULL,"-verbose",&verbose));
 32:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test BVGetSplitRows (length %" PetscInt_FMT ", l=%" PetscInt_FMT ", k=%" PetscInt_FMT ", nc=%" PetscInt_FMT ").\n",n,l,k,nc));

 34:   /* Create Nest matrix */
 35:   PetscCall(MatCreate(PETSC_COMM_WORLD,&T));
 36:   PetscCall(MatSetSizes(T,PETSC_DECIDE,PETSC_DECIDE,n,n));
 37:   PetscCall(MatSetFromOptions(T));
 38:   PetscCall(MatGetOwnershipRange(T,&rstart,&rend));
 39:   for (i=rstart;i<rend;i++) {
 40:     if (i>0) PetscCall(MatSetValue(T,i,i-1,-1.0,INSERT_VALUES));
 41:     PetscCall(MatSetValue(T,i,i,2.0,INSERT_VALUES));
 42:     if (i<n-1) PetscCall(MatSetValue(T,i,i+1,-1.0,INSERT_VALUES));
 43:   }
 44:   PetscCall(MatAssemblyBegin(T,MAT_FINAL_ASSEMBLY));
 45:   PetscCall(MatAssemblyEnd(T,MAT_FINAL_ASSEMBLY));
 46:   PetscCall(MatCreateConstantDiagonal(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,n,n,1.0,&D));

 48:   block[0] = T;
 49:   block[1] = block[2] = NULL;
 50:   block[3] = D;
 51:   PetscCall(MatCreateNest(PETSC_COMM_WORLD,2,NULL,2,NULL,block,&M));
 52:   PetscCall(MatDestroy(&T));
 53:   PetscCall(MatDestroy(&D));
 54:   PetscCall(MatNestGetISs(M,is,NULL));

 56:   PetscCall(MatView(M,NULL));

 58:   /* Create template vector */
 59:   PetscCall(MatCreateVecs(M,&t,NULL));

 61:   /* Create BV object X */
 62:   PetscCall(BVCreate(PETSC_COMM_WORLD,&X));
 63:   PetscCall(PetscObjectSetName((PetscObject)X,"X"));
 64:   PetscCall(BVSetSizesFromVec(X,t,k));
 65:   PetscCall(BVSetFromOptions(X));

 67:   /* Generate constraints and attach them to X */
 68:   if (nc>0) {
 69:     PetscCall(VecDuplicateVecs(t,nc,&C));
 70:     for (j=0;j<nc;j++) {
 71:       for (i=0;i<=j;i++) PetscCall(VecSetValue(C[j],nc-i+1,1.0,INSERT_VALUES));
 72:       PetscCall(VecAssemblyBegin(C[j]));
 73:       PetscCall(VecAssemblyEnd(C[j]));
 74:     }
 75:     PetscCall(BVInsertConstraints(X,&nc,C));
 76:     PetscCall(VecDestroyVecs(nc,&C));
 77:   }

 79:   /* Set up viewer */
 80:   PetscCall(PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&view));
 81:   if (verbose) PetscCall(PetscViewerPushFormat(view,PETSC_VIEWER_ASCII_MATLAB));

 83:   /* Fill X entries */
 84:   for (j=0;j<k;j++) {
 85:     PetscCall(BVGetColumn(X,j,&v));
 86:     PetscCall(VecSet(v,0.0));
 87:     PetscCall(VecGetSubVector(v,is[0],&v1));
 88:     PetscCall(VecGetSubVector(v,is[1],&v2));
 89:     for (i=0;i<4;i++) {
 90:       if (i+j>=rstart && i+j<rend) {
 91:         PetscCall(VecSetValue(v1,i+j,(PetscScalar)(3*i+j-2),INSERT_VALUES));
 92:         PetscCall(VecSetValue(v2,i+j,(PetscScalar)(-i+2*(j+1)),INSERT_VALUES));
 93:       }
 94:     }
 95:     PetscCall(VecAssemblyBegin(v1));
 96:     PetscCall(VecAssemblyBegin(v2));
 97:     PetscCall(VecAssemblyEnd(v1));
 98:     PetscCall(VecAssemblyEnd(v2));
 99:     PetscCall(VecRestoreSubVector(v,is[0],&v1));
100:     PetscCall(VecRestoreSubVector(v,is[1],&v2));
101:     PetscCall(BVRestoreColumn(X,j,&v));
102:   }
103:   if (verbose) PetscCall(BVView(X,view));

105:   /* Get split BVs */
106:   PetscCall(BVGetSplitRows(X,is[0],is[1],&U,&L));
107:   PetscCall(PetscObjectSetName((PetscObject)U,"U"));
108:   PetscCall(PetscObjectSetName((PetscObject)L,"L"));

110:   if (verbose) {
111:     PetscCall(BVView(U,view));
112:     PetscCall(BVView(L,view));
113:   }

115:   /* Copy l-th column of U to first column of L */
116:   PetscCall(BVGetColumn(U,l,&v));
117:   PetscCall(BVGetColumn(L,0,&w));
118:   PetscCall(VecCopy(v,w));
119:   PetscCall(BVRestoreColumn(U,l,&v));
120:   PetscCall(BVRestoreColumn(L,0,&w));

122:   /* Finished using the split BVs */
123:   PetscCall(BVRestoreSplitRows(X,is[0],is[1],&U,&L));
124:   if (verbose) PetscCall(BVView(X,view));

126:   /* Check: print bottom part of first column */
127:   PetscCall(BVGetColumn(X,0,&v));
128:   PetscCall(VecGetSubVector(v,is[1],&v2));
129:   PetscCall(VecView(v2,NULL));
130:   PetscCall(VecRestoreSubVector(v,is[1],&v2));
131:   PetscCall(BVRestoreColumn(X,0,&v));

133:   PetscCall(BVDestroy(&X));
134:   PetscCall(VecDestroy(&t));
135:   PetscCall(MatDestroy(&M));
136:   PetscCall(SlepcFinalize());
137:   return 0;
138: }

140: /*TEST

142:    testset:
143:       nsize: {{1 2}}
144:       output_file: output/test19_1.out
145:       filter: grep -v Process | grep -v Object | sed -e 's/mpi/seq/' | sed -e 's/seqcuda/seq/' | sed -e 's/seqaijcusparse/seqaij/' | sed -e 's/seqhip/seq/' | sed -e 's/seqaijhipsparse/seqaij/' | sed -e 's/nc=2/nc=0/'
146:       test:
147:          suffix: 1
148:          args: -nc {{0 2}} -bv_type {{svec mat}}
149:       test:
150:          suffix: 1_cuda
151:          args: -nc {{0 2}} -bv_type {{svec mat}} -mat_type aijcusparse
152:          requires: cuda
153:       test:
154:          suffix: 1_hip
155:          args: -nc {{0 2}} -bv_type {{svec mat}} -mat_type aijhipsparse
156:          requires: hip

158: TEST*/