Actual source code: submatfree.c

  1: #include <petsctao.h>
  2: #include <../src/tao/matrix/submatfree.h>

  4: /*@C
  5:   MatCreateSubMatrixFree - Creates a reduced matrix by masking a
  6:   full matrix.

  8:    Collective on matrix

 10:    Input Parameters:
 11: +  mat - matrix of arbitrary type
 12: .  Rows - the rows that will be in the submatrix
 13: -  Cols - the columns that will be in the submatrix

 15:    Output Parameters:
 16: .  J - New matrix

 18:    Notes:
 19:    The caller is responsible for destroying the input objects after matrix J has been destroyed.

 21:    Level: developer

 23: .seealso: MatCreate()
 24: @*/
 25: PetscErrorCode MatCreateSubMatrixFree(Mat mat,IS Rows, IS Cols, Mat *J)
 26: {
 27:   MPI_Comm         comm=PetscObjectComm((PetscObject)mat);
 28:   MatSubMatFreeCtx ctx;
 29:   PetscErrorCode   ierr;
 30:   PetscInt         mloc,nloc,m,n;

 33:   PetscNew(&ctx);
 34:   ctx->A  = mat;
 35:   MatGetSize(mat,&m,&n);
 36:   MatGetLocalSize(mat,&mloc,&nloc);
 37:   MatCreateVecs(mat,NULL,&ctx->VC);
 38:   ctx->VR = ctx->VC;
 39:    PetscObjectReference((PetscObject)mat);

 41:   ctx->Rows = Rows;
 42:   ctx->Cols = Cols;
 43:   PetscObjectReference((PetscObject)Rows);
 44:   PetscObjectReference((PetscObject)Cols);
 45:   MatCreateShell(comm,mloc,nloc,m,n,ctx,J);
 46:   MatShellSetManageScalingShifts(*J);
 47:   MatShellSetOperation(*J,MATOP_MULT,(void(*)(void))MatMult_SMF);
 48:   MatShellSetOperation(*J,MATOP_DESTROY,(void(*)(void))MatDestroy_SMF);
 49:   MatShellSetOperation(*J,MATOP_VIEW,(void(*)(void))MatView_SMF);
 50:   MatShellSetOperation(*J,MATOP_MULT_TRANSPOSE,(void(*)(void))MatMultTranspose_SMF);
 51:   MatShellSetOperation(*J,MATOP_DIAGONAL_SET,(void(*)(void))MatDiagonalSet_SMF);
 52:   MatShellSetOperation(*J,MATOP_SHIFT,(void(*)(void))MatShift_SMF);
 53:   MatShellSetOperation(*J,MATOP_EQUAL,(void(*)(void))MatEqual_SMF);
 54:   MatShellSetOperation(*J,MATOP_SCALE,(void(*)(void))MatScale_SMF);
 55:   MatShellSetOperation(*J,MATOP_TRANSPOSE,(void(*)(void))MatTranspose_SMF);
 56:   MatShellSetOperation(*J,MATOP_GET_DIAGONAL,(void(*)(void))MatGetDiagonal_SMF);
 57:   MatShellSetOperation(*J,MATOP_CREATE_SUBMATRICES,(void(*)(void))MatCreateSubMatrices_SMF);
 58:   MatShellSetOperation(*J,MATOP_NORM,(void(*)(void))MatNorm_SMF);
 59:   MatShellSetOperation(*J,MATOP_DUPLICATE,(void(*)(void))MatDuplicate_SMF);
 60:   MatShellSetOperation(*J,MATOP_CREATE_SUBMATRIX,(void(*)(void))MatCreateSubMatrix_SMF);
 61:   MatShellSetOperation(*J,MATOP_GET_ROW_MAX,(void(*)(void))MatDuplicate_SMF);

 63:   PetscLogObjectParent((PetscObject)mat,(PetscObject)(*J));
 64:   return(0);
 65: }

 67: PetscErrorCode MatSMFResetRowColumn(Mat mat,IS Rows,IS Cols)
 68: {
 69:   MatSubMatFreeCtx ctx;
 70:   PetscErrorCode   ierr;

 73:   MatShellGetContext(mat,&ctx);
 74:   ISDestroy(&ctx->Rows);
 75:   ISDestroy(&ctx->Cols);
 76:   PetscObjectReference((PetscObject)Rows);
 77:   PetscObjectReference((PetscObject)Cols);
 78:   ctx->Cols=Cols;
 79:   ctx->Rows=Rows;
 80:   return(0);
 81: }

 83: PetscErrorCode MatMult_SMF(Mat mat,Vec a,Vec y)
 84: {
 85:   MatSubMatFreeCtx ctx;
 86:   PetscErrorCode   ierr;

 89:   MatShellGetContext(mat,&ctx);
 90:   VecCopy(a,ctx->VR);
 91:   VecISSet(ctx->VR,ctx->Cols,0.0);
 92:   MatMult(ctx->A,ctx->VR,y);
 93:   VecISSet(y,ctx->Rows,0.0);
 94:   return(0);
 95: }

 97: PetscErrorCode MatMultTranspose_SMF(Mat mat,Vec a,Vec y)
 98: {
 99:   MatSubMatFreeCtx ctx;
100:   PetscErrorCode   ierr;

103:   MatShellGetContext(mat,&ctx);
104:   VecCopy(a,ctx->VC);
105:   VecISSet(ctx->VC,ctx->Rows,0.0);
106:   MatMultTranspose(ctx->A,ctx->VC,y);
107:   VecISSet(y,ctx->Cols,0.0);
108:   return(0);
109: }

111: PetscErrorCode MatDiagonalSet_SMF(Mat M, Vec D,InsertMode is)
112: {
113:   MatSubMatFreeCtx ctx;
114:   PetscErrorCode   ierr;

117:   MatShellGetContext(M,&ctx);
118:   MatDiagonalSet(ctx->A,D,is);
119:   return(0);
120: }

122: PetscErrorCode MatDestroy_SMF(Mat mat)
123: {
124:   PetscErrorCode   ierr;
125:   MatSubMatFreeCtx ctx;

128:   MatShellGetContext(mat,&ctx);
129:   MatDestroy(&ctx->A);
130:   ISDestroy(&ctx->Rows);
131:   ISDestroy(&ctx->Cols);
132:   VecDestroy(&ctx->VC);
133:   PetscFree(ctx);
134:   return(0);
135: }

137: PetscErrorCode MatView_SMF(Mat mat,PetscViewer viewer)
138: {
139:   PetscErrorCode   ierr;
140:   MatSubMatFreeCtx ctx;

143:   MatShellGetContext(mat,&ctx);
144:   MatView(ctx->A,viewer);
145:   return(0);
146: }

148: PetscErrorCode MatShift_SMF(Mat Y, PetscReal a)
149: {
150:   PetscErrorCode   ierr;
151:   MatSubMatFreeCtx ctx;

154:   MatShellGetContext(Y,&ctx);
155:   MatShift(ctx->A,a);
156:   return(0);
157: }

159: PetscErrorCode MatDuplicate_SMF(Mat mat,MatDuplicateOption op,Mat *M)
160: {
161:   PetscErrorCode   ierr;
162:   MatSubMatFreeCtx ctx;

165:   MatShellGetContext(mat,&ctx);
166:   MatCreateSubMatrixFree(ctx->A,ctx->Rows,ctx->Cols,M);
167:   return(0);
168: }

170: PetscErrorCode MatEqual_SMF(Mat A,Mat B,PetscBool *flg)
171: {
172:   PetscErrorCode    ierr;
173:   MatSubMatFreeCtx  ctx1,ctx2;
174:   PetscBool         flg1,flg2,flg3;

177:   MatShellGetContext(A,&ctx1);
178:   MatShellGetContext(B,&ctx2);
179:   ISEqual(ctx1->Rows,ctx2->Rows,&flg2);
180:   ISEqual(ctx1->Cols,ctx2->Cols,&flg3);
181:   if (flg2==PETSC_FALSE || flg3==PETSC_FALSE) {
182:     *flg=PETSC_FALSE;
183:   } else {
184:     MatEqual(ctx1->A,ctx2->A,&flg1);
185:     if (flg1==PETSC_FALSE) { *flg=PETSC_FALSE;}
186:     else { *flg=PETSC_TRUE;}
187:   }
188:   return(0);
189: }

191: PetscErrorCode MatScale_SMF(Mat mat, PetscReal a)
192: {
193:   PetscErrorCode   ierr;
194:   MatSubMatFreeCtx ctx;

197:   MatShellGetContext(mat,&ctx);
198:   MatScale(ctx->A,a);
199:   return(0);
200: }

202: PetscErrorCode MatTranspose_SMF(Mat mat,Mat *B)
203: {
205:   PetscFunctionReturn(1);
206: }

208: PetscErrorCode MatGetDiagonal_SMF(Mat mat,Vec v)
209: {
210:   PetscErrorCode   ierr;
211:   MatSubMatFreeCtx ctx;

214:   MatShellGetContext(mat,&ctx);
215:   MatGetDiagonal(ctx->A,v);
216:   return(0);
217: }

219: PetscErrorCode MatGetRowMax_SMF(Mat M, Vec D)
220: {
221:   MatSubMatFreeCtx ctx;
222:   PetscErrorCode   ierr;

225:   MatShellGetContext(M,&ctx);
226:   MatGetRowMax(ctx->A,D,NULL);
227:   return(0);
228: }

230: PetscErrorCode MatCreateSubMatrices_SMF(Mat A,PetscInt n, IS *irow,IS *icol,MatReuse scall,Mat **B)
231: {
233:   PetscInt       i;

236:   if (scall == MAT_INITIAL_MATRIX) {
237:     PetscCalloc1(n+1,B);
238:   }

240:   for (i=0; i<n; i++) {
241:     MatCreateSubMatrix_SMF(A,irow[i],icol[i],scall,&(*B)[i]);
242:   }
243:   return(0);
244: }

246: PetscErrorCode MatCreateSubMatrix_SMF(Mat mat,IS isrow,IS iscol,MatReuse cll,
247:                         Mat *newmat)
248: {
249:   PetscErrorCode   ierr;
250:   MatSubMatFreeCtx ctx;

253:   MatShellGetContext(mat,&ctx);
254:   if (newmat) {
255:     MatDestroy(&*newmat);
256:   }
257:   MatCreateSubMatrixFree(ctx->A,isrow,iscol, newmat);
258:   return(0);
259: }

261: PetscErrorCode MatGetRow_SMF(Mat mat,PetscInt row,PetscInt *ncols,const PetscInt **cols,const PetscScalar **vals)
262: {
263:   PetscErrorCode   ierr;
264:   MatSubMatFreeCtx ctx;

267:   MatShellGetContext(mat,&ctx);
268:   MatGetRow(ctx->A,row,ncols,cols,vals);
269:   return(0);
270: }

272: PetscErrorCode MatRestoreRow_SMF(Mat mat,PetscInt row,PetscInt *ncols,const PetscInt **cols,const PetscScalar **vals)
273: {
274:   PetscErrorCode   ierr;
275:   MatSubMatFreeCtx ctx;

278:   MatShellGetContext(mat,&ctx);
279:   MatRestoreRow(ctx->A,row,ncols,cols,vals);
280:   return(0);
281: }

283: PetscErrorCode MatGetColumnVector_SMF(Mat mat,Vec Y, PetscInt col)
284: {
285:   PetscErrorCode   ierr;
286:   MatSubMatFreeCtx ctx;

289:   MatShellGetContext(mat,&ctx);
290:   MatGetColumnVector(ctx->A,Y,col);
291:   return(0);
292: }

294: PetscErrorCode MatNorm_SMF(Mat mat,NormType type,PetscReal *norm)
295: {
296:   PetscErrorCode    ierr;
297:   MatSubMatFreeCtx  ctx;

300:   MatShellGetContext(mat,&ctx);
301:   if (type == NORM_FROBENIUS) {
302:     *norm = 1.0;
303:   } else if (type == NORM_1 || type == NORM_INFINITY) {
304:     *norm = 1.0;
305:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No two norm");
306:   return(0);
307: }