Actual source code: ex97.c

petsc-3.9.4 2018-09-11
Report Typos and Errors
  1: static const char help[] = "Tests MatCreateSubMatrix with MatSubMatrix versus MatAIJ, non-square\n";

  3:  #include <petscmat.h>

  5: static PetscErrorCode AssembleMatrix(MPI_Comm comm,Mat *A)
  6: {
  8:   Mat            B;
  9:   PetscInt       i,ms,me;

 12:   MatCreate(comm,&B);
 13:   MatSetSizes(B,5,6,PETSC_DETERMINE,PETSC_DETERMINE);
 14:   MatSetFromOptions(B);
 15:   MatSetUp(B);
 16:   MatGetOwnershipRange(B,&ms,&me);
 17:   for (i=ms; i<me; i++) {
 18:     MatSetValue(B,i,i,1.0*i,INSERT_VALUES);
 19:   }
 20:   MatSetValue(B,me-1,me,me*me,INSERT_VALUES);
 21:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
 22:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
 23:   *A   = B;
 24:   return(0);
 25: }

 27: static PetscErrorCode Compare2(Vec *X,const char *test)
 28: {
 30:   PetscReal      norm;
 31:   Vec            Y;
 32:   PetscInt       verbose = 0;

 35:   VecDuplicate(X[0],&Y);
 36:   VecCopy(X[0],Y);
 37:   VecAYPX(Y,-1.0,X[1]);
 38:   VecNorm(Y,NORM_INFINITY,&norm);

 40:   PetscOptionsGetInt(NULL,NULL,"-verbose",&verbose,NULL);
 41:   if (norm < 1.e-12 && verbose < 1) {
 42:     PetscPrintf(PETSC_COMM_WORLD,"%30s: norm difference < 1e-12\n",test);
 43:   } else {
 44:     PetscPrintf(PETSC_COMM_WORLD,"%30s: norm difference %g\n",test,(double)norm);
 45:   }
 46:   if (verbose > 1) {
 47:     VecView(X[0],PETSC_VIEWER_STDOUT_WORLD);
 48:     VecView(X[1],PETSC_VIEWER_STDOUT_WORLD);
 49:     VecView(Y,PETSC_VIEWER_STDOUT_WORLD);
 50:   }
 51:   VecDestroy(&Y);
 52:   return(0);
 53: }

 55: static PetscErrorCode CheckMatrices(Mat A,Mat B,Vec left,Vec right,Vec X,Vec Y,Vec X1,Vec Y1)
 56: {
 58:   Vec            *ltmp,*rtmp;

 61:   VecDuplicateVecs(right,2,&rtmp);
 62:   VecDuplicateVecs(left,2,&ltmp);
 63:   MatScale(A,PETSC_PI);
 64:   MatScale(B,PETSC_PI);
 65:   MatDiagonalScale(A,left,right);
 66:   MatDiagonalScale(B,left,right);

 68:   MatMult(A,X,ltmp[0]);
 69:   MatMult(B,X,ltmp[1]);
 70:   Compare2(ltmp,"MatMult");

 72:   MatMultTranspose(A,Y,rtmp[0]);
 73:   MatMultTranspose(B,Y,rtmp[1]);
 74:   Compare2(rtmp,"MatMultTranspose");

 76:   VecCopy(Y1,ltmp[0]);
 77:   VecCopy(Y1,ltmp[1]);
 78:   MatMultAdd(A,X,ltmp[0],ltmp[0]);
 79:   MatMultAdd(B,X,ltmp[1],ltmp[1]);
 80:   Compare2(ltmp,"MatMultAdd v2==v3");

 82:   MatMultAdd(A,X,Y1,ltmp[0]);
 83:   MatMultAdd(B,X,Y1,ltmp[1]);
 84:   Compare2(ltmp,"MatMultAdd v2!=v3");

 86:   VecCopy(X1,rtmp[0]);
 87:   VecCopy(X1,rtmp[1]);
 88:   MatMultTransposeAdd(A,Y,rtmp[0],rtmp[0]);
 89:   MatMultTransposeAdd(B,Y,rtmp[1],rtmp[1]);
 90:   Compare2(rtmp,"MatMultTransposeAdd v2==v3");

 92:   MatMultTransposeAdd(A,Y,X1,rtmp[0]);
 93:   MatMultTransposeAdd(B,Y,X1,rtmp[1]);
 94:   Compare2(rtmp,"MatMultTransposeAdd v2!=v3");

 96:   VecDestroyVecs(2,&ltmp);
 97:   VecDestroyVecs(2,&rtmp);
 98:   return(0);
 99: }

101: int main(int argc, char *argv[])
102: {
104:   Mat            A,B,Asub,Bsub;
105:   PetscInt       ms,idxrow[3],idxcol[4];
106:   Vec            left,right,X,Y,X1,Y1;
107:   IS             isrow,iscol;
108:   PetscBool      random = PETSC_TRUE;

110:   PetscInitialize(&argc,&argv,NULL,help);
111:   AssembleMatrix(PETSC_COMM_WORLD,&A);
112:   AssembleMatrix(PETSC_COMM_WORLD,&B);
113:   MatSetOperation(B,MATOP_CREATE_SUBMATRIX,NULL);
114:   MatSetOperation(B,MATOP_CREATE_SUBMATRICES,NULL);
115:   MatGetOwnershipRange(A,&ms,NULL);

117:   idxrow[0] = ms+1;
118:   idxrow[1] = ms+2;
119:   idxrow[2] = ms+4;
120:   ISCreateGeneral(PETSC_COMM_WORLD,3,idxrow,PETSC_USE_POINTER,&isrow);

122:   idxcol[0] = ms+1;
123:   idxcol[1] = ms+2;
124:   idxcol[2] = ms+4;
125:   idxcol[3] = ms+5;
126:   ISCreateGeneral(PETSC_COMM_WORLD,4,idxcol,PETSC_USE_POINTER,&iscol);

128:   MatCreateSubMatrix(A,isrow,iscol,MAT_INITIAL_MATRIX,&Asub);
129:   MatCreateSubMatrix(B,isrow,iscol,MAT_INITIAL_MATRIX,&Bsub);

131:   MatCreateVecs(Asub,&right,&left);
132:   VecDuplicate(right,&X);
133:   VecDuplicate(right,&X1);
134:   VecDuplicate(left,&Y);
135:   VecDuplicate(left,&Y1);

137:   PetscOptionsGetBool(NULL,NULL,"-random",&random,NULL);
138:   if (random) {
139:     VecSetRandom(right,NULL);
140:     VecSetRandom(left,NULL);
141:     VecSetRandom(X,NULL);
142:     VecSetRandom(Y,NULL);
143:     VecSetRandom(X1,NULL);
144:     VecSetRandom(Y1,NULL);
145:   } else {
146:     VecSet(right,1.0);
147:     VecSet(left,2.0);
148:     VecSet(X,3.0);
149:     VecSet(Y,4.0);
150:     VecSet(X1,3.0);
151:     VecSet(Y1,4.0);
152:   }
153:   CheckMatrices(Asub,Bsub,left,right,X,Y,X1,Y1);
154:   ISDestroy(&isrow);
155:   ISDestroy(&iscol);
156:   MatDestroy(&A);
157:   MatDestroy(&B);
158:   MatDestroy(&Asub);
159:   MatDestroy(&Bsub);
160:   VecDestroy(&left);
161:   VecDestroy(&right);
162:   VecDestroy(&X);
163:   VecDestroy(&Y);
164:   VecDestroy(&X1);
165:   VecDestroy(&Y1);
166:   PetscFinalize();
167:   return ierr;
168: }



172: /*TEST

174:    test:
175:       nsize: 3

177: TEST*/