1 #include <petsctao.h> /*I "petsctao.h" I*/ 2 #include <../src/tao/matrix/submatfree.h> /*I "submatfree.h" I*/ 3 4 /*@C 5 MatCreateSubMatrixFree - Creates a reduced matrix by masking a 6 full matrix. 7 8 Collective on matrix 9 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 14 15 Output Parameters: 16 . J - New matrix 17 18 Notes: 19 The caller is responsible for destroying the input objects after matrix J has been destroyed. 20 21 Level: developer 22 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; 31 32 PetscFunctionBegin; 33 ierr = PetscNew(&ctx);CHKERRQ(ierr); 34 ctx->A = mat; 35 ierr = MatGetSize(mat,&m,&n);CHKERRQ(ierr); 36 ierr = MatGetLocalSize(mat,&mloc,&nloc);CHKERRQ(ierr); 37 ierr = MatCreateVecs(mat,NULL,&ctx->VC);CHKERRQ(ierr); 38 ctx->VR = ctx->VC; 39 ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 40 41 ctx->Rows = Rows; 42 ctx->Cols = Cols; 43 ierr = PetscObjectReference((PetscObject)Rows);CHKERRQ(ierr); 44 ierr = PetscObjectReference((PetscObject)Cols);CHKERRQ(ierr); 45 ierr = MatCreateShell(comm,mloc,nloc,m,n,ctx,J);CHKERRQ(ierr); 46 ierr = MatShellSetManageScalingShifts(*J);CHKERRQ(ierr); 47 ierr = MatShellSetOperation(*J,MATOP_MULT,(void(*)(void))MatMult_SMF);CHKERRQ(ierr); 48 ierr = MatShellSetOperation(*J,MATOP_DESTROY,(void(*)(void))MatDestroy_SMF);CHKERRQ(ierr); 49 ierr = MatShellSetOperation(*J,MATOP_VIEW,(void(*)(void))MatView_SMF);CHKERRQ(ierr); 50 ierr = MatShellSetOperation(*J,MATOP_MULT_TRANSPOSE,(void(*)(void))MatMultTranspose_SMF);CHKERRQ(ierr); 51 ierr = MatShellSetOperation(*J,MATOP_DIAGONAL_SET,(void(*)(void))MatDiagonalSet_SMF);CHKERRQ(ierr); 52 ierr = MatShellSetOperation(*J,MATOP_SHIFT,(void(*)(void))MatShift_SMF);CHKERRQ(ierr); 53 ierr = MatShellSetOperation(*J,MATOP_EQUAL,(void(*)(void))MatEqual_SMF);CHKERRQ(ierr); 54 ierr = MatShellSetOperation(*J,MATOP_SCALE,(void(*)(void))MatScale_SMF);CHKERRQ(ierr); 55 ierr = MatShellSetOperation(*J,MATOP_TRANSPOSE,(void(*)(void))MatTranspose_SMF);CHKERRQ(ierr); 56 ierr = MatShellSetOperation(*J,MATOP_GET_DIAGONAL,(void(*)(void))MatGetDiagonal_SMF);CHKERRQ(ierr); 57 ierr = MatShellSetOperation(*J,MATOP_CREATE_SUBMATRICES,(void(*)(void))MatCreateSubMatrices_SMF);CHKERRQ(ierr); 58 ierr = MatShellSetOperation(*J,MATOP_NORM,(void(*)(void))MatNorm_SMF);CHKERRQ(ierr); 59 ierr = MatShellSetOperation(*J,MATOP_DUPLICATE,(void(*)(void))MatDuplicate_SMF);CHKERRQ(ierr); 60 ierr = MatShellSetOperation(*J,MATOP_CREATE_SUBMATRIX,(void(*)(void))MatCreateSubMatrix_SMF);CHKERRQ(ierr); 61 ierr = MatShellSetOperation(*J,MATOP_GET_ROW_MAX,(void(*)(void))MatDuplicate_SMF);CHKERRQ(ierr); 62 63 ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)(*J));CHKERRQ(ierr); 64 PetscFunctionReturn(0); 65 } 66 67 PetscErrorCode MatSMFResetRowColumn(Mat mat,IS Rows,IS Cols) 68 { 69 MatSubMatFreeCtx ctx; 70 PetscErrorCode ierr; 71 72 PetscFunctionBegin; 73 ierr = MatShellGetContext(mat,&ctx);CHKERRQ(ierr); 74 ierr = ISDestroy(&ctx->Rows);CHKERRQ(ierr); 75 ierr = ISDestroy(&ctx->Cols);CHKERRQ(ierr); 76 ierr = PetscObjectReference((PetscObject)Rows);CHKERRQ(ierr); 77 ierr = PetscObjectReference((PetscObject)Cols);CHKERRQ(ierr); 78 ctx->Cols=Cols; 79 ctx->Rows=Rows; 80 PetscFunctionReturn(0); 81 } 82 83 PetscErrorCode MatMult_SMF(Mat mat,Vec a,Vec y) 84 { 85 MatSubMatFreeCtx ctx; 86 PetscErrorCode ierr; 87 88 PetscFunctionBegin; 89 ierr = MatShellGetContext(mat,&ctx);CHKERRQ(ierr); 90 ierr = VecCopy(a,ctx->VR);CHKERRQ(ierr); 91 ierr = VecISSet(ctx->VR,ctx->Cols,0.0);CHKERRQ(ierr); 92 ierr = MatMult(ctx->A,ctx->VR,y);CHKERRQ(ierr); 93 ierr = VecISSet(y,ctx->Rows,0.0);CHKERRQ(ierr); 94 PetscFunctionReturn(0); 95 } 96 97 PetscErrorCode MatMultTranspose_SMF(Mat mat,Vec a,Vec y) 98 { 99 MatSubMatFreeCtx ctx; 100 PetscErrorCode ierr; 101 102 PetscFunctionBegin; 103 ierr = MatShellGetContext(mat,&ctx);CHKERRQ(ierr); 104 ierr = VecCopy(a,ctx->VC);CHKERRQ(ierr); 105 ierr = VecISSet(ctx->VC,ctx->Rows,0.0);CHKERRQ(ierr); 106 ierr = MatMultTranspose(ctx->A,ctx->VC,y);CHKERRQ(ierr); 107 ierr = VecISSet(y,ctx->Cols,0.0);CHKERRQ(ierr); 108 PetscFunctionReturn(0); 109 } 110 111 PetscErrorCode MatDiagonalSet_SMF(Mat M, Vec D,InsertMode is) 112 { 113 MatSubMatFreeCtx ctx; 114 PetscErrorCode ierr; 115 116 PetscFunctionBegin; 117 ierr = MatShellGetContext(M,&ctx);CHKERRQ(ierr); 118 ierr = MatDiagonalSet(ctx->A,D,is);CHKERRQ(ierr); 119 PetscFunctionReturn(0); 120 } 121 122 PetscErrorCode MatDestroy_SMF(Mat mat) 123 { 124 PetscErrorCode ierr; 125 MatSubMatFreeCtx ctx; 126 127 PetscFunctionBegin; 128 ierr = MatShellGetContext(mat,&ctx);CHKERRQ(ierr); 129 ierr = MatDestroy(&ctx->A);CHKERRQ(ierr); 130 ierr = ISDestroy(&ctx->Rows);CHKERRQ(ierr); 131 ierr = ISDestroy(&ctx->Cols);CHKERRQ(ierr); 132 ierr = VecDestroy(&ctx->VC);CHKERRQ(ierr); 133 ierr = PetscFree(ctx);CHKERRQ(ierr); 134 PetscFunctionReturn(0); 135 } 136 137 PetscErrorCode MatView_SMF(Mat mat,PetscViewer viewer) 138 { 139 PetscErrorCode ierr; 140 MatSubMatFreeCtx ctx; 141 142 PetscFunctionBegin; 143 ierr = MatShellGetContext(mat,&ctx);CHKERRQ(ierr); 144 ierr = MatView(ctx->A,viewer);CHKERRQ(ierr); 145 PetscFunctionReturn(0); 146 } 147 148 PetscErrorCode MatShift_SMF(Mat Y, PetscReal a) 149 { 150 PetscErrorCode ierr; 151 MatSubMatFreeCtx ctx; 152 153 PetscFunctionBegin; 154 ierr = MatShellGetContext(Y,&ctx);CHKERRQ(ierr); 155 ierr = MatShift(ctx->A,a);CHKERRQ(ierr); 156 PetscFunctionReturn(0); 157 } 158 159 PetscErrorCode MatDuplicate_SMF(Mat mat,MatDuplicateOption op,Mat *M) 160 { 161 PetscErrorCode ierr; 162 MatSubMatFreeCtx ctx; 163 164 PetscFunctionBegin; 165 ierr = MatShellGetContext(mat,&ctx);CHKERRQ(ierr); 166 ierr = MatCreateSubMatrixFree(ctx->A,ctx->Rows,ctx->Cols,M);CHKERRQ(ierr); 167 PetscFunctionReturn(0); 168 } 169 170 PetscErrorCode MatEqual_SMF(Mat A,Mat B,PetscBool *flg) 171 { 172 PetscErrorCode ierr; 173 MatSubMatFreeCtx ctx1,ctx2; 174 PetscBool flg1,flg2,flg3; 175 176 PetscFunctionBegin; 177 ierr = MatShellGetContext(A,&ctx1);CHKERRQ(ierr); 178 ierr = MatShellGetContext(B,&ctx2);CHKERRQ(ierr); 179 ierr = ISEqual(ctx1->Rows,ctx2->Rows,&flg2);CHKERRQ(ierr); 180 ierr = ISEqual(ctx1->Cols,ctx2->Cols,&flg3);CHKERRQ(ierr); 181 if (flg2==PETSC_FALSE || flg3==PETSC_FALSE) { 182 *flg=PETSC_FALSE; 183 } else { 184 ierr = MatEqual(ctx1->A,ctx2->A,&flg1);CHKERRQ(ierr); 185 if (flg1==PETSC_FALSE) { *flg=PETSC_FALSE;} 186 else { *flg=PETSC_TRUE;} 187 } 188 PetscFunctionReturn(0); 189 } 190 191 PetscErrorCode MatScale_SMF(Mat mat, PetscReal a) 192 { 193 PetscErrorCode ierr; 194 MatSubMatFreeCtx ctx; 195 196 PetscFunctionBegin; 197 ierr = MatShellGetContext(mat,&ctx);CHKERRQ(ierr); 198 ierr = MatScale(ctx->A,a);CHKERRQ(ierr); 199 PetscFunctionReturn(0); 200 } 201 202 PetscErrorCode MatTranspose_SMF(Mat mat,Mat *B) 203 { 204 PetscFunctionBegin; 205 PetscFunctionReturn(1); 206 } 207 208 PetscErrorCode MatGetDiagonal_SMF(Mat mat,Vec v) 209 { 210 PetscErrorCode ierr; 211 MatSubMatFreeCtx ctx; 212 213 PetscFunctionBegin; 214 ierr = MatShellGetContext(mat,&ctx);CHKERRQ(ierr); 215 ierr = MatGetDiagonal(ctx->A,v);CHKERRQ(ierr); 216 PetscFunctionReturn(0); 217 } 218 219 PetscErrorCode MatGetRowMax_SMF(Mat M, Vec D) 220 { 221 MatSubMatFreeCtx ctx; 222 PetscErrorCode ierr; 223 224 PetscFunctionBegin; 225 ierr = MatShellGetContext(M,&ctx);CHKERRQ(ierr); 226 ierr = MatGetRowMax(ctx->A,D,NULL);CHKERRQ(ierr); 227 PetscFunctionReturn(0); 228 } 229 230 PetscErrorCode MatCreateSubMatrices_SMF(Mat A,PetscInt n, IS *irow,IS *icol,MatReuse scall,Mat **B) 231 { 232 PetscErrorCode ierr; 233 PetscInt i; 234 235 PetscFunctionBegin; 236 if (scall == MAT_INITIAL_MATRIX) { 237 ierr = PetscCalloc1(n+1,B);CHKERRQ(ierr); 238 } 239 240 for (i=0; i<n; i++) { 241 ierr = MatCreateSubMatrix_SMF(A,irow[i],icol[i],scall,&(*B)[i]);CHKERRQ(ierr); 242 } 243 PetscFunctionReturn(0); 244 } 245 246 PetscErrorCode MatCreateSubMatrix_SMF(Mat mat,IS isrow,IS iscol,MatReuse cll, 247 Mat *newmat) 248 { 249 PetscErrorCode ierr; 250 MatSubMatFreeCtx ctx; 251 252 PetscFunctionBegin; 253 ierr = MatShellGetContext(mat,&ctx);CHKERRQ(ierr); 254 if (newmat) { 255 ierr = MatDestroy(&*newmat);CHKERRQ(ierr); 256 } 257 ierr = MatCreateSubMatrixFree(ctx->A,isrow,iscol, newmat);CHKERRQ(ierr); 258 PetscFunctionReturn(0); 259 } 260 261 PetscErrorCode MatGetRow_SMF(Mat mat,PetscInt row,PetscInt *ncols,const PetscInt **cols,const PetscScalar **vals) 262 { 263 PetscErrorCode ierr; 264 MatSubMatFreeCtx ctx; 265 266 PetscFunctionBegin; 267 ierr = MatShellGetContext(mat,&ctx);CHKERRQ(ierr); 268 ierr = MatGetRow(ctx->A,row,ncols,cols,vals);CHKERRQ(ierr); 269 PetscFunctionReturn(0); 270 } 271 272 PetscErrorCode MatRestoreRow_SMF(Mat mat,PetscInt row,PetscInt *ncols,const PetscInt **cols,const PetscScalar **vals) 273 { 274 PetscErrorCode ierr; 275 MatSubMatFreeCtx ctx; 276 277 PetscFunctionBegin; 278 ierr = MatShellGetContext(mat,&ctx);CHKERRQ(ierr); 279 ierr = MatRestoreRow(ctx->A,row,ncols,cols,vals);CHKERRQ(ierr); 280 PetscFunctionReturn(0); 281 } 282 283 PetscErrorCode MatGetColumnVector_SMF(Mat mat,Vec Y, PetscInt col) 284 { 285 PetscErrorCode ierr; 286 MatSubMatFreeCtx ctx; 287 288 PetscFunctionBegin; 289 ierr = MatShellGetContext(mat,&ctx);CHKERRQ(ierr); 290 ierr = MatGetColumnVector(ctx->A,Y,col);CHKERRQ(ierr); 291 PetscFunctionReturn(0); 292 } 293 294 PetscErrorCode MatNorm_SMF(Mat mat,NormType type,PetscReal *norm) 295 { 296 PetscErrorCode ierr; 297 MatSubMatFreeCtx ctx; 298 299 PetscFunctionBegin; 300 ierr = MatShellGetContext(mat,&ctx);CHKERRQ(ierr); 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 PetscFunctionReturn(0); 307 } 308 309