1be1d678aSKris Buschelman 2b4319ba4SBarry Smith /* 3b4319ba4SBarry Smith Creates a matrix class for using the Neumann-Neumann type preconditioners. 4b4319ba4SBarry Smith This stores the matrices in globally unassembled form. Each processor 5b4319ba4SBarry Smith assembles only its local Neumann problem and the parallel matrix vector 6b4319ba4SBarry Smith product is handled "implicitly". 7b4319ba4SBarry Smith 8b4319ba4SBarry Smith We provide: 9b4319ba4SBarry Smith MatMult() 10b4319ba4SBarry Smith 11b4319ba4SBarry Smith Currently this allows for only one subdomain per processor. 12b4319ba4SBarry Smith 13b4319ba4SBarry Smith */ 14b4319ba4SBarry Smith 15c6db04a5SJed Brown #include <../src/mat/impls/is/matis.h> /*I "petscmat.h" I*/ 16b4319ba4SBarry Smith 17b4319ba4SBarry Smith #undef __FUNCT__ 18*4e4c7dbeSStefano Zampini #define __FUNCT__ "MatGetVecs_IS" 19*4e4c7dbeSStefano Zampini PetscErrorCode MatGetVecs_IS(Mat A,Vec *right,Vec *left) 20*4e4c7dbeSStefano Zampini { 21*4e4c7dbeSStefano Zampini PetscErrorCode ierr; 22*4e4c7dbeSStefano Zampini PetscInt bs,rows,columns; 23*4e4c7dbeSStefano Zampini 24*4e4c7dbeSStefano Zampini PetscFunctionBegin; 25*4e4c7dbeSStefano Zampini ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr); 26*4e4c7dbeSStefano Zampini ierr = MatGetSize(A,&rows,&columns);CHKERRQ(ierr); 27*4e4c7dbeSStefano Zampini if(right) { 28*4e4c7dbeSStefano Zampini ierr = VecCreate(((PetscObject)A)->comm,right);CHKERRQ(ierr); 29*4e4c7dbeSStefano Zampini ierr = VecSetBlockSize(*right,bs);CHKERRQ(ierr); 30*4e4c7dbeSStefano Zampini ierr = VecSetSizes(*right,PETSC_DECIDE,columns);CHKERRQ(ierr); 31*4e4c7dbeSStefano Zampini ierr = VecSetType(*right,VECMPI);CHKERRQ(ierr); 32*4e4c7dbeSStefano Zampini } 33*4e4c7dbeSStefano Zampini if(left) { 34*4e4c7dbeSStefano Zampini ierr = VecCreate(((PetscObject)A)->comm,left);CHKERRQ(ierr); 35*4e4c7dbeSStefano Zampini ierr = VecSetBlockSize(*left,bs);CHKERRQ(ierr); 36*4e4c7dbeSStefano Zampini ierr = VecSetSizes(*left,PETSC_DECIDE,rows);CHKERRQ(ierr); 37*4e4c7dbeSStefano Zampini ierr = VecSetType(*left,VECMPI);CHKERRQ(ierr); 38*4e4c7dbeSStefano Zampini } 39*4e4c7dbeSStefano Zampini PetscFunctionReturn(0); 40*4e4c7dbeSStefano Zampini } 41*4e4c7dbeSStefano Zampini 42*4e4c7dbeSStefano Zampini #undef __FUNCT__ 43b4319ba4SBarry Smith #define __FUNCT__ "MatDestroy_IS" 44dfbe8321SBarry Smith PetscErrorCode MatDestroy_IS(Mat A) 45b4319ba4SBarry Smith { 46dfbe8321SBarry Smith PetscErrorCode ierr; 47b4319ba4SBarry Smith Mat_IS *b = (Mat_IS*)A->data; 48b4319ba4SBarry Smith 49b4319ba4SBarry Smith PetscFunctionBegin; 506bf464f9SBarry Smith ierr = MatDestroy(&b->A);CHKERRQ(ierr); 516bf464f9SBarry Smith ierr = VecScatterDestroy(&b->ctx);CHKERRQ(ierr); 526bf464f9SBarry Smith ierr = VecDestroy(&b->x);CHKERRQ(ierr); 536bf464f9SBarry Smith ierr = VecDestroy(&b->y);CHKERRQ(ierr); 546bf464f9SBarry Smith ierr = ISLocalToGlobalMappingDestroy(&b->mapping);CHKERRQ(ierr); 55bf0cc555SLisandro Dalcin ierr = PetscFree(A->data);CHKERRQ(ierr); 56dbd8c25aSHong Zhang ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr); 57901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C","",PETSC_NULL);CHKERRQ(ierr); 58b4319ba4SBarry Smith PetscFunctionReturn(0); 59b4319ba4SBarry Smith } 60b4319ba4SBarry Smith 61b4319ba4SBarry Smith #undef __FUNCT__ 62b4319ba4SBarry Smith #define __FUNCT__ "MatMult_IS" 63dfbe8321SBarry Smith PetscErrorCode MatMult_IS(Mat A,Vec x,Vec y) 64b4319ba4SBarry Smith { 65dfbe8321SBarry Smith PetscErrorCode ierr; 66b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 67b4319ba4SBarry Smith PetscScalar zero = 0.0; 68b4319ba4SBarry Smith 69b4319ba4SBarry Smith PetscFunctionBegin; 70b4319ba4SBarry Smith /* scatter the global vector x into the local work vector */ 71ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 72ca9f406cSSatish Balay ierr = VecScatterEnd(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 73b4319ba4SBarry Smith 74b4319ba4SBarry Smith /* multiply the local matrix */ 75b4319ba4SBarry Smith ierr = MatMult(is->A,is->x,is->y);CHKERRQ(ierr); 76b4319ba4SBarry Smith 77b4319ba4SBarry Smith /* scatter product back into global memory */ 782dcb1b2aSMatthew Knepley ierr = VecSet(y,zero);CHKERRQ(ierr); 79ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 80ca9f406cSSatish Balay ierr = VecScatterEnd(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 81b4319ba4SBarry Smith 82b4319ba4SBarry Smith PetscFunctionReturn(0); 83b4319ba4SBarry Smith } 84b4319ba4SBarry Smith 85b4319ba4SBarry Smith #undef __FUNCT__ 862e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultAdd_IS" 872e74eeadSLisandro Dalcin static PetscErrorCode MatMultAdd_IS(Mat A,Vec v1,Vec v2,Vec v3) 882e74eeadSLisandro Dalcin { 892e74eeadSLisandro Dalcin PetscErrorCode ierr; 902e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 912e74eeadSLisandro Dalcin 922e74eeadSLisandro Dalcin PetscFunctionBegin; /* v3 = v2 + A * v1.*/ 932e74eeadSLisandro Dalcin /* scatter the global vector v1 into the local work vector */ 94ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,v1,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 95ca9f406cSSatish Balay ierr = VecScatterEnd (is->ctx,v1,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 96ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,v2,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 97ca9f406cSSatish Balay ierr = VecScatterEnd (is->ctx,v2,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 982e74eeadSLisandro Dalcin 992e74eeadSLisandro Dalcin /* multiply the local matrix */ 1002e74eeadSLisandro Dalcin ierr = MatMultAdd(is->A,is->x,is->y,is->y);CHKERRQ(ierr); 1012e74eeadSLisandro Dalcin 1022e74eeadSLisandro Dalcin /* scatter result back into global vector */ 1032e74eeadSLisandro Dalcin ierr = VecSet(v3,0);CHKERRQ(ierr); 104ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,is->y,v3,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 105ca9f406cSSatish Balay ierr = VecScatterEnd (is->ctx,is->y,v3,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1062e74eeadSLisandro Dalcin PetscFunctionReturn(0); 1072e74eeadSLisandro Dalcin } 1082e74eeadSLisandro Dalcin 1092e74eeadSLisandro Dalcin #undef __FUNCT__ 1102e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTranspose_IS" 1112e74eeadSLisandro Dalcin PetscErrorCode MatMultTranspose_IS(Mat A,Vec x,Vec y) 1122e74eeadSLisandro Dalcin { 1132e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 1142e74eeadSLisandro Dalcin PetscErrorCode ierr; 1152e74eeadSLisandro Dalcin 1162e74eeadSLisandro Dalcin PetscFunctionBegin; /* y = A' * x */ 1172e74eeadSLisandro Dalcin /* scatter the global vector x into the local work vector */ 118ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 119ca9f406cSSatish Balay ierr = VecScatterEnd(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1202e74eeadSLisandro Dalcin 1212e74eeadSLisandro Dalcin /* multiply the local matrix */ 1222e74eeadSLisandro Dalcin ierr = MatMultTranspose(is->A,is->x,is->y);CHKERRQ(ierr); 1232e74eeadSLisandro Dalcin 1242e74eeadSLisandro Dalcin /* scatter product back into global vector */ 1252e74eeadSLisandro Dalcin ierr = VecSet(y,0);CHKERRQ(ierr); 126ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 127ca9f406cSSatish Balay ierr = VecScatterEnd(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1282e74eeadSLisandro Dalcin PetscFunctionReturn(0); 1292e74eeadSLisandro Dalcin } 1302e74eeadSLisandro Dalcin 1312e74eeadSLisandro Dalcin #undef __FUNCT__ 1322e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTransposeAdd_IS" 1332e74eeadSLisandro Dalcin PetscErrorCode MatMultTransposeAdd_IS(Mat A,Vec v1,Vec v2,Vec v3) 1342e74eeadSLisandro Dalcin { 1352e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 1362e74eeadSLisandro Dalcin PetscErrorCode ierr; 1372e74eeadSLisandro Dalcin 1382e74eeadSLisandro Dalcin PetscFunctionBegin; /* v3 = v2 + A' * v1.*/ 1392e74eeadSLisandro Dalcin /* scatter the global vectors v1 and v2 into the local work vectors */ 140ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,v1,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 141ca9f406cSSatish Balay ierr = VecScatterEnd (is->ctx,v1,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 142ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,v2,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 143ca9f406cSSatish Balay ierr = VecScatterEnd (is->ctx,v2,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1442e74eeadSLisandro Dalcin 1452e74eeadSLisandro Dalcin /* multiply the local matrix */ 1462e74eeadSLisandro Dalcin ierr = MatMultTransposeAdd(is->A,is->x,is->y,is->y);CHKERRQ(ierr); 1472e74eeadSLisandro Dalcin 1482e74eeadSLisandro Dalcin /* scatter result back into global vector */ 1492e74eeadSLisandro Dalcin ierr = VecSet(v3,0);CHKERRQ(ierr); 150ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,is->y,v3,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 151ca9f406cSSatish Balay ierr = VecScatterEnd (is->ctx,is->y,v3,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1522e74eeadSLisandro Dalcin PetscFunctionReturn(0); 1532e74eeadSLisandro Dalcin } 1542e74eeadSLisandro Dalcin 1552e74eeadSLisandro Dalcin #undef __FUNCT__ 156b4319ba4SBarry Smith #define __FUNCT__ "MatView_IS" 157dfbe8321SBarry Smith PetscErrorCode MatView_IS(Mat A,PetscViewer viewer) 158b4319ba4SBarry Smith { 159b4319ba4SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 160dfbe8321SBarry Smith PetscErrorCode ierr; 161b4319ba4SBarry Smith PetscViewer sviewer; 162b4319ba4SBarry Smith 163b4319ba4SBarry Smith PetscFunctionBegin; 164b4319ba4SBarry Smith ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr); 165b4319ba4SBarry Smith ierr = MatView(a->A,sviewer);CHKERRQ(ierr); 166b4319ba4SBarry Smith ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr); 167b4319ba4SBarry Smith PetscFunctionReturn(0); 168b4319ba4SBarry Smith } 169b4319ba4SBarry Smith 170b4319ba4SBarry Smith #undef __FUNCT__ 171b4319ba4SBarry Smith #define __FUNCT__ "MatSetLocalToGlobalMapping_IS" 172784ac674SJed Brown PetscErrorCode MatSetLocalToGlobalMapping_IS(Mat A,ISLocalToGlobalMapping rmapping,ISLocalToGlobalMapping cmapping) 173b4319ba4SBarry Smith { 174dfbe8321SBarry Smith PetscErrorCode ierr; 175*4e4c7dbeSStefano Zampini PetscInt n,bs; 176b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 177b4319ba4SBarry Smith IS from,to; 178b4319ba4SBarry Smith Vec global; 179b4319ba4SBarry Smith 180b4319ba4SBarry Smith PetscFunctionBegin; 181e7e72b3dSBarry Smith if (is->mapping) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Mapping already set for matrix"); 182784ac674SJed Brown PetscCheckSameComm(A,1,rmapping,2); 183784ac674SJed Brown if (rmapping != cmapping) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_INCOMP,"MATIS requires the row and column mappings to be identical"); 184784ac674SJed Brown ierr = PetscObjectReference((PetscObject)rmapping);CHKERRQ(ierr); 1856bf464f9SBarry Smith ierr = ISLocalToGlobalMappingDestroy(&is->mapping);CHKERRQ(ierr); 186784ac674SJed Brown is->mapping = rmapping; 187b4319ba4SBarry Smith 188b4319ba4SBarry Smith /* Create the local matrix A */ 189784ac674SJed Brown ierr = ISLocalToGlobalMappingGetSize(rmapping,&n);CHKERRQ(ierr); 190*4e4c7dbeSStefano Zampini ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr); 191f69a0ea3SMatthew Knepley ierr = MatCreate(PETSC_COMM_SELF,&is->A);CHKERRQ(ierr); 192f69a0ea3SMatthew Knepley ierr = MatSetSizes(is->A,n,n,n,n);CHKERRQ(ierr); 193*4e4c7dbeSStefano Zampini ierr = MatSetBlockSize(is->A,bs);CHKERRQ(ierr); 1942e74eeadSLisandro Dalcin ierr = MatSetOptionsPrefix(is->A,"is");CHKERRQ(ierr); 195b4319ba4SBarry Smith ierr = MatSetFromOptions(is->A);CHKERRQ(ierr); 196b4319ba4SBarry Smith 197b4319ba4SBarry Smith /* Create the local work vectors */ 198*4e4c7dbeSStefano Zampini ierr = VecCreate(PETSC_COMM_SELF,&is->x);CHKERRQ(ierr); 199*4e4c7dbeSStefano Zampini ierr = VecSetBlockSize(is->x,bs);CHKERRQ(ierr); 200*4e4c7dbeSStefano Zampini ierr = VecSetSizes(is->x,n,n);CHKERRQ(ierr); 201*4e4c7dbeSStefano Zampini ierr = VecSetOptionsPrefix(is->x,"is");CHKERRQ(ierr); 202*4e4c7dbeSStefano Zampini ierr = VecSetFromOptions(is->x);CHKERRQ(ierr); 203b4319ba4SBarry Smith ierr = VecDuplicate(is->x,&is->y);CHKERRQ(ierr); 204b4319ba4SBarry Smith 205b4319ba4SBarry Smith /* setup the global to local scatter */ 206b4319ba4SBarry Smith ierr = ISCreateStride(PETSC_COMM_SELF,n,0,1,&to);CHKERRQ(ierr); 207784ac674SJed Brown ierr = ISLocalToGlobalMappingApplyIS(rmapping,to,&from);CHKERRQ(ierr); 208*4e4c7dbeSStefano Zampini ierr = MatGetVecs(A,&global,PETSC_NULL);CHKERRQ(ierr); 209b4319ba4SBarry Smith ierr = VecScatterCreate(global,from,is->x,to,&is->ctx);CHKERRQ(ierr); 2106bf464f9SBarry Smith ierr = VecDestroy(&global);CHKERRQ(ierr); 2116bf464f9SBarry Smith ierr = ISDestroy(&to);CHKERRQ(ierr); 2126bf464f9SBarry Smith ierr = ISDestroy(&from);CHKERRQ(ierr); 213b4319ba4SBarry Smith PetscFunctionReturn(0); 214b4319ba4SBarry Smith } 215b4319ba4SBarry Smith 2162e74eeadSLisandro Dalcin #define ISG2LMapApply(mapping,n,in,out) 0;\ 2172e74eeadSLisandro Dalcin if (!(mapping)->globals) {\ 2182e74eeadSLisandro Dalcin PetscErrorCode _ierr = ISGlobalToLocalMappingApply((mapping),IS_GTOLM_MASK,0,0,0,0);CHKERRQ(_ierr);\ 2192e74eeadSLisandro Dalcin }\ 2202e74eeadSLisandro Dalcin {\ 2212e74eeadSLisandro Dalcin PetscInt _i,*_globals = (mapping)->globals,_start = (mapping)->globalstart,_end = (mapping)->globalend;\ 2222e74eeadSLisandro Dalcin for (_i=0; _i<n; _i++) {\ 2232e74eeadSLisandro Dalcin if (in[_i] < 0) out[_i] = in[_i];\ 2242e74eeadSLisandro Dalcin else if (in[_i] < _start) out[_i] = -1;\ 2252e74eeadSLisandro Dalcin else if (in[_i] > _end) out[_i] = -1;\ 2262e74eeadSLisandro Dalcin else out[_i] = _globals[in[_i] - _start];\ 2272e74eeadSLisandro Dalcin }\ 2282e74eeadSLisandro Dalcin } 2292e74eeadSLisandro Dalcin 2302e74eeadSLisandro Dalcin #undef __FUNCT__ 2312e74eeadSLisandro Dalcin #define __FUNCT__ "MatSetValues_IS" 2322e74eeadSLisandro Dalcin PetscErrorCode MatSetValues_IS(Mat mat, PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols, const PetscScalar *values, InsertMode addv) 2332e74eeadSLisandro Dalcin { 2342e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)mat->data; 2352e74eeadSLisandro Dalcin PetscInt rows_l[2048],cols_l[2048]; 2362e74eeadSLisandro Dalcin PetscErrorCode ierr; 2372e74eeadSLisandro Dalcin 2382e74eeadSLisandro Dalcin PetscFunctionBegin; 2392e74eeadSLisandro Dalcin #if defined(PETSC_USE_DEBUG) 2402e74eeadSLisandro Dalcin if (m > 2048 || n > 2048) { 241e32f2f54SBarry Smith SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Number of row/column indices must be <= 2048: they are %D %D",m,n); 2422e74eeadSLisandro Dalcin } 2432e74eeadSLisandro Dalcin #endif 2442e74eeadSLisandro Dalcin ierr = ISG2LMapApply(is->mapping,m,rows,rows_l);CHKERRQ(ierr); 2452e74eeadSLisandro Dalcin ierr = ISG2LMapApply(is->mapping,n,cols,cols_l);CHKERRQ(ierr); 2462e74eeadSLisandro Dalcin ierr = MatSetValues(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr); 2472e74eeadSLisandro Dalcin PetscFunctionReturn(0); 2482e74eeadSLisandro Dalcin } 2492e74eeadSLisandro Dalcin 2502e74eeadSLisandro Dalcin #undef ISG2LMapSetUp 2512e74eeadSLisandro Dalcin #undef ISG2LMapApply 252b4319ba4SBarry Smith 253b4319ba4SBarry Smith #undef __FUNCT__ 254b4319ba4SBarry Smith #define __FUNCT__ "MatSetValuesLocal_IS" 25513f74950SBarry Smith PetscErrorCode MatSetValuesLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv) 256b4319ba4SBarry Smith { 257dfbe8321SBarry Smith PetscErrorCode ierr; 258b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 259b4319ba4SBarry Smith 260b4319ba4SBarry Smith PetscFunctionBegin; 261b4319ba4SBarry Smith ierr = MatSetValues(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr); 262b4319ba4SBarry Smith PetscFunctionReturn(0); 263b4319ba4SBarry Smith } 264b4319ba4SBarry Smith 265b4319ba4SBarry Smith #undef __FUNCT__ 2662e74eeadSLisandro Dalcin #define __FUNCT__ "MatZeroRows_IS" 2672b40b63fSBarry Smith PetscErrorCode MatZeroRows_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b) 2682e74eeadSLisandro Dalcin { 2692e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 2702e74eeadSLisandro Dalcin PetscInt n_l=0, *rows_l = PETSC_NULL; 2712e74eeadSLisandro Dalcin PetscErrorCode ierr; 2722e74eeadSLisandro Dalcin 2732e74eeadSLisandro Dalcin PetscFunctionBegin; 27497b48c8fSBarry Smith if (x && b) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"No support"); 2752e74eeadSLisandro Dalcin if (n) { 2762e74eeadSLisandro Dalcin ierr = PetscMalloc(n*sizeof(PetscInt),&rows_l);CHKERRQ(ierr); 2772e74eeadSLisandro Dalcin ierr = ISGlobalToLocalMappingApply(is->mapping,IS_GTOLM_DROP,n,rows,&n_l,rows_l);CHKERRQ(ierr); 2782e74eeadSLisandro Dalcin } 2792b40b63fSBarry Smith ierr = MatZeroRowsLocal(A,n_l,rows_l,diag,x,b);CHKERRQ(ierr); 280c31cb41cSBarry Smith ierr = PetscFree(rows_l);CHKERRQ(ierr); 2812e74eeadSLisandro Dalcin PetscFunctionReturn(0); 2822e74eeadSLisandro Dalcin } 2832e74eeadSLisandro Dalcin 2842e74eeadSLisandro Dalcin #undef __FUNCT__ 285b4319ba4SBarry Smith #define __FUNCT__ "MatZeroRowsLocal_IS" 2862b40b63fSBarry Smith PetscErrorCode MatZeroRowsLocal_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b) 287b4319ba4SBarry Smith { 288b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 289dfbe8321SBarry Smith PetscErrorCode ierr; 290f4df32b1SMatthew Knepley PetscInt i; 291b4319ba4SBarry Smith PetscScalar *array; 292b4319ba4SBarry Smith 293b4319ba4SBarry Smith PetscFunctionBegin; 2949c3e2b52SBarry Smith if (x && b) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"No support"); 295b4319ba4SBarry Smith { 296b4319ba4SBarry Smith /* 297b4319ba4SBarry Smith Set up is->x as a "counting vector". This is in order to MatMult_IS 298b4319ba4SBarry Smith work properly in the interface nodes. 299b4319ba4SBarry Smith */ 300b4319ba4SBarry Smith Vec counter; 301b4319ba4SBarry Smith PetscScalar one=1.0, zero=0.0; 302*4e4c7dbeSStefano Zampini ierr = MatGetVecs(A,&counter,PETSC_NULL);CHKERRQ(ierr); 3032dcb1b2aSMatthew Knepley ierr = VecSet(counter,zero);CHKERRQ(ierr); 3042dcb1b2aSMatthew Knepley ierr = VecSet(is->x,one);CHKERRQ(ierr); 305ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,is->x,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 306ca9f406cSSatish Balay ierr = VecScatterEnd (is->ctx,is->x,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 307ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,counter,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 308ca9f406cSSatish Balay ierr = VecScatterEnd (is->ctx,counter,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 3096bf464f9SBarry Smith ierr = VecDestroy(&counter);CHKERRQ(ierr); 310b4319ba4SBarry Smith } 311958c9bccSBarry Smith if (!n) { 312b4319ba4SBarry Smith is->pure_neumann = PETSC_TRUE; 313b4319ba4SBarry Smith } else { 314b4319ba4SBarry Smith is->pure_neumann = PETSC_FALSE; 315b4319ba4SBarry Smith ierr = VecGetArray(is->x,&array);CHKERRQ(ierr); 3162b40b63fSBarry Smith ierr = MatZeroRows(is->A,n,rows,diag,0,0);CHKERRQ(ierr); 317b4319ba4SBarry Smith for (i=0; i<n; i++) { 318f4df32b1SMatthew Knepley ierr = MatSetValue(is->A,rows[i],rows[i],diag/(array[rows[i]]),INSERT_VALUES);CHKERRQ(ierr); 319b4319ba4SBarry Smith } 320b4319ba4SBarry Smith ierr = MatAssemblyBegin(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 321b4319ba4SBarry Smith ierr = MatAssemblyEnd (is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 322b4319ba4SBarry Smith ierr = VecRestoreArray(is->x,&array);CHKERRQ(ierr); 323b4319ba4SBarry Smith } 324b4319ba4SBarry Smith 325b4319ba4SBarry Smith PetscFunctionReturn(0); 326b4319ba4SBarry Smith } 327b4319ba4SBarry Smith 328b4319ba4SBarry Smith #undef __FUNCT__ 329b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyBegin_IS" 330dfbe8321SBarry Smith PetscErrorCode MatAssemblyBegin_IS(Mat A,MatAssemblyType type) 331b4319ba4SBarry Smith { 332b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 333dfbe8321SBarry Smith PetscErrorCode ierr; 334b4319ba4SBarry Smith 335b4319ba4SBarry Smith PetscFunctionBegin; 336b4319ba4SBarry Smith ierr = MatAssemblyBegin(is->A,type);CHKERRQ(ierr); 337b4319ba4SBarry Smith PetscFunctionReturn(0); 338b4319ba4SBarry Smith } 339b4319ba4SBarry Smith 340b4319ba4SBarry Smith #undef __FUNCT__ 341b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyEnd_IS" 342dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_IS(Mat A,MatAssemblyType type) 343b4319ba4SBarry Smith { 344b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 345dfbe8321SBarry Smith PetscErrorCode ierr; 346b4319ba4SBarry Smith 347b4319ba4SBarry Smith PetscFunctionBegin; 348b4319ba4SBarry Smith ierr = MatAssemblyEnd(is->A,type);CHKERRQ(ierr); 349b4319ba4SBarry Smith PetscFunctionReturn(0); 350b4319ba4SBarry Smith } 351b4319ba4SBarry Smith 352b4319ba4SBarry Smith EXTERN_C_BEGIN 353b4319ba4SBarry Smith #undef __FUNCT__ 354b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat_IS" 3557087cfbeSBarry Smith PetscErrorCode MatISGetLocalMat_IS(Mat mat,Mat *local) 356b4319ba4SBarry Smith { 357b4319ba4SBarry Smith Mat_IS *is = (Mat_IS *)mat->data; 358b4319ba4SBarry Smith 359b4319ba4SBarry Smith PetscFunctionBegin; 360b4319ba4SBarry Smith *local = is->A; 361b4319ba4SBarry Smith PetscFunctionReturn(0); 362b4319ba4SBarry Smith } 363b4319ba4SBarry Smith EXTERN_C_END 364b4319ba4SBarry Smith 365b4319ba4SBarry Smith #undef __FUNCT__ 366b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat" 367b4319ba4SBarry Smith /*@ 368b4319ba4SBarry Smith MatISGetLocalMat - Gets the local matrix stored inside a MATIS matrix. 369b4319ba4SBarry Smith 370b4319ba4SBarry Smith Input Parameter: 371b4319ba4SBarry Smith . mat - the matrix 372b4319ba4SBarry Smith 373b4319ba4SBarry Smith Output Parameter: 374b4319ba4SBarry Smith . local - the local matrix usually MATSEQAIJ 375b4319ba4SBarry Smith 376b4319ba4SBarry Smith Level: advanced 377b4319ba4SBarry Smith 378b4319ba4SBarry Smith Notes: 379b4319ba4SBarry Smith This can be called if you have precomputed the nonzero structure of the 380b4319ba4SBarry Smith matrix and want to provide it to the inner matrix object to improve the performance 381b4319ba4SBarry Smith of the MatSetValues() operation. 382b4319ba4SBarry Smith 383b4319ba4SBarry Smith .seealso: MATIS 384b4319ba4SBarry Smith @*/ 3857087cfbeSBarry Smith PetscErrorCode MatISGetLocalMat(Mat mat,Mat *local) 386b4319ba4SBarry Smith { 3874ac538c5SBarry Smith PetscErrorCode ierr; 388b4319ba4SBarry Smith 389b4319ba4SBarry Smith PetscFunctionBegin; 3900700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 391b4319ba4SBarry Smith PetscValidPointer(local,2); 3924ac538c5SBarry Smith ierr = PetscUseMethod(mat,"MatISGetLocalMat_C",(Mat,Mat *),(mat,local));CHKERRQ(ierr); 393b4319ba4SBarry Smith PetscFunctionReturn(0); 394b4319ba4SBarry Smith } 395b4319ba4SBarry Smith 3963b03a366Sstefano_zampini EXTERN_C_BEGIN 3973b03a366Sstefano_zampini #undef __FUNCT__ 3983b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat_IS" 3993b03a366Sstefano_zampini PetscErrorCode MatISSetLocalMat_IS(Mat mat,Mat local) 4003b03a366Sstefano_zampini { 4013b03a366Sstefano_zampini Mat_IS *is = (Mat_IS *)mat->data; 4023b03a366Sstefano_zampini PetscInt nrows,ncols,orows,ocols; 4033b03a366Sstefano_zampini PetscErrorCode ierr; 4043b03a366Sstefano_zampini 4053b03a366Sstefano_zampini PetscFunctionBegin; 406*4e4c7dbeSStefano Zampini if(is->A) { 4073b03a366Sstefano_zampini ierr = MatGetSize(is->A,&orows,&ocols);CHKERRQ(ierr); 4083b03a366Sstefano_zampini ierr = MatGetSize(local,&nrows,&ncols);CHKERRQ(ierr); 4093b03a366Sstefano_zampini if(orows != nrows || ocols != ncols ) { 4103b03a366Sstefano_zampini SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Local MATIS matrix should be of size %dx%d (you passed a %dx%d matrix)\n",orows,ocols,nrows,ncols); 4113b03a366Sstefano_zampini } 412*4e4c7dbeSStefano Zampini } 4133b03a366Sstefano_zampini ierr = PetscObjectReference((PetscObject)local);CHKERRQ(ierr); 4143b03a366Sstefano_zampini ierr = MatDestroy(&is->A);CHKERRQ(ierr); 4153b03a366Sstefano_zampini is->A = local; 4163b03a366Sstefano_zampini 4173b03a366Sstefano_zampini PetscFunctionReturn(0); 4183b03a366Sstefano_zampini } 4193b03a366Sstefano_zampini EXTERN_C_END 4203b03a366Sstefano_zampini 4213b03a366Sstefano_zampini #undef __FUNCT__ 4223b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat" 4233b03a366Sstefano_zampini /*@ 4243b03a366Sstefano_zampini MatISSetLocalMat - Set the local matrix stored inside a MATIS matrix. 4253b03a366Sstefano_zampini 4263b03a366Sstefano_zampini Input Parameter: 4273b03a366Sstefano_zampini . mat - the matrix 4283b03a366Sstefano_zampini . local - the local matrix usually MATSEQAIJ 4293b03a366Sstefano_zampini 4303b03a366Sstefano_zampini Output Parameter: 4313b03a366Sstefano_zampini 4323b03a366Sstefano_zampini Level: advanced 4333b03a366Sstefano_zampini 4343b03a366Sstefano_zampini Notes: 4353b03a366Sstefano_zampini This can be called if you have precomputed the local matrix and 4363b03a366Sstefano_zampini want to provide it to the matrix object MATIS. 4373b03a366Sstefano_zampini 4383b03a366Sstefano_zampini .seealso: MATIS 4393b03a366Sstefano_zampini @*/ 4403b03a366Sstefano_zampini PetscErrorCode MatISSetLocalMat(Mat mat,Mat local) 4413b03a366Sstefano_zampini { 4423b03a366Sstefano_zampini PetscErrorCode ierr; 4433b03a366Sstefano_zampini 4443b03a366Sstefano_zampini PetscFunctionBegin; 4453b03a366Sstefano_zampini PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 4463b03a366Sstefano_zampini ierr = PetscUseMethod(mat,"MatISSetLocalMat_C",(Mat,Mat),(mat,local));CHKERRQ(ierr); 4473b03a366Sstefano_zampini PetscFunctionReturn(0); 4483b03a366Sstefano_zampini } 4493b03a366Sstefano_zampini 4506726f965SBarry Smith #undef __FUNCT__ 4516726f965SBarry Smith #define __FUNCT__ "MatZeroEntries_IS" 4526726f965SBarry Smith PetscErrorCode MatZeroEntries_IS(Mat A) 4536726f965SBarry Smith { 4546726f965SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 4556726f965SBarry Smith PetscErrorCode ierr; 4566726f965SBarry Smith 4576726f965SBarry Smith PetscFunctionBegin; 4586726f965SBarry Smith ierr = MatZeroEntries(a->A);CHKERRQ(ierr); 4596726f965SBarry Smith PetscFunctionReturn(0); 4606726f965SBarry Smith } 4616726f965SBarry Smith 4626726f965SBarry Smith #undef __FUNCT__ 4632e74eeadSLisandro Dalcin #define __FUNCT__ "MatScale_IS" 4642e74eeadSLisandro Dalcin PetscErrorCode MatScale_IS(Mat A,PetscScalar a) 4652e74eeadSLisandro Dalcin { 4662e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 4672e74eeadSLisandro Dalcin PetscErrorCode ierr; 4682e74eeadSLisandro Dalcin 4692e74eeadSLisandro Dalcin PetscFunctionBegin; 4702e74eeadSLisandro Dalcin ierr = MatScale(is->A,a);CHKERRQ(ierr); 4712e74eeadSLisandro Dalcin PetscFunctionReturn(0); 4722e74eeadSLisandro Dalcin } 4732e74eeadSLisandro Dalcin 4742e74eeadSLisandro Dalcin #undef __FUNCT__ 4752e74eeadSLisandro Dalcin #define __FUNCT__ "MatGetDiagonal_IS" 4762e74eeadSLisandro Dalcin PetscErrorCode MatGetDiagonal_IS(Mat A, Vec v) 4772e74eeadSLisandro Dalcin { 4782e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 4792e74eeadSLisandro Dalcin PetscErrorCode ierr; 4802e74eeadSLisandro Dalcin 4812e74eeadSLisandro Dalcin PetscFunctionBegin; 4822e74eeadSLisandro Dalcin /* get diagonal of the local matrix */ 4832e74eeadSLisandro Dalcin ierr = MatGetDiagonal(is->A,is->x);CHKERRQ(ierr); 4842e74eeadSLisandro Dalcin 4852e74eeadSLisandro Dalcin /* scatter diagonal back into global vector */ 4862e74eeadSLisandro Dalcin ierr = VecSet(v,0);CHKERRQ(ierr); 487ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,is->x,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 488ca9f406cSSatish Balay ierr = VecScatterEnd (is->ctx,is->x,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 4892e74eeadSLisandro Dalcin PetscFunctionReturn(0); 4902e74eeadSLisandro Dalcin } 4912e74eeadSLisandro Dalcin 4922e74eeadSLisandro Dalcin #undef __FUNCT__ 4936726f965SBarry Smith #define __FUNCT__ "MatSetOption_IS" 494ace3abfcSBarry Smith PetscErrorCode MatSetOption_IS(Mat A,MatOption op,PetscBool flg) 4956726f965SBarry Smith { 4966726f965SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 4976726f965SBarry Smith PetscErrorCode ierr; 4986726f965SBarry Smith 4996726f965SBarry Smith PetscFunctionBegin; 5004e0d8c25SBarry Smith ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr); 5016726f965SBarry Smith PetscFunctionReturn(0); 5026726f965SBarry Smith } 5036726f965SBarry Smith 504284134d9SBarry Smith #undef __FUNCT__ 505284134d9SBarry Smith #define __FUNCT__ "MatCreateIS" 506284134d9SBarry Smith /*@ 507284134d9SBarry Smith MatCreateIS - Creates a "process" unassmembled matrix, it is assembled on each 508284134d9SBarry Smith process but not across processes. 509284134d9SBarry Smith 510284134d9SBarry Smith Input Parameters: 511284134d9SBarry Smith + comm - MPI communicator that will share the matrix 512*4e4c7dbeSStefano Zampini . bs - local and global block size of the matrix 513284134d9SBarry Smith . m,n,M,N - local and/or global sizes of the the left and right vector used in matrix vector products 514284134d9SBarry Smith - map - mapping that defines the global number for each local number 515284134d9SBarry Smith 516284134d9SBarry Smith Output Parameter: 517284134d9SBarry Smith . A - the resulting matrix 518284134d9SBarry Smith 5198e6c10adSSatish Balay Level: advanced 5208e6c10adSSatish Balay 521284134d9SBarry Smith Notes: See MATIS for more details 5228cda6cd7SBarry Smith m and n are NOT related to the size of the map, they are the size of the part of the vector owned 5238cda6cd7SBarry Smith by that process. m + nghosts (or n + nghosts) is the length of map since map maps all local points 5248cda6cd7SBarry Smith plus the ghost points to global indices. 525284134d9SBarry Smith 526284134d9SBarry Smith .seealso: MATIS, MatSetLocalToGlobalMapping() 527284134d9SBarry Smith @*/ 528*4e4c7dbeSStefano Zampini PetscErrorCode MatCreateIS(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,ISLocalToGlobalMapping map,Mat *A) 529284134d9SBarry Smith { 530284134d9SBarry Smith PetscErrorCode ierr; 531*4e4c7dbeSStefano Zampini PetscInt local_size; 532*4e4c7dbeSStefano Zampini Vec global; 533284134d9SBarry Smith 534284134d9SBarry Smith PetscFunctionBegin; 535284134d9SBarry Smith ierr = MatCreate(comm,A);CHKERRQ(ierr); 536284134d9SBarry Smith ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr); 537284134d9SBarry Smith ierr = MatSetType(*A,MATIS);CHKERRQ(ierr); 538*4e4c7dbeSStefano Zampini /* Set block size manually */ 539*4e4c7dbeSStefano Zampini (*A)->rmap->bs=bs; 540*4e4c7dbeSStefano Zampini (*A)->cmap->bs=bs; 5413b03a366Sstefano_zampini ierr = MatSetUp(*A);CHKERRQ(ierr); 542784ac674SJed Brown ierr = MatSetLocalToGlobalMapping(*A,map,map);CHKERRQ(ierr); 543*4e4c7dbeSStefano Zampini /* force local dimensions */ 544*4e4c7dbeSStefano Zampini ierr = MatGetVecs(*A,&global,PETSC_NULL);CHKERRQ(ierr); 545*4e4c7dbeSStefano Zampini ierr = VecGetLocalSize(global,&local_size);CHKERRQ(ierr); 546*4e4c7dbeSStefano Zampini (*A)->rmap->n=local_size; 547*4e4c7dbeSStefano Zampini (*A)->cmap->n=local_size; 548284134d9SBarry Smith PetscFunctionReturn(0); 549284134d9SBarry Smith } 550284134d9SBarry Smith 551b4319ba4SBarry Smith /*MC 552b4319ba4SBarry Smith MATIS - MATIS = "is" - A matrix type to be used for using the Neumann-Neumann type preconditioners. 553b4319ba4SBarry Smith This stores the matrices in globally unassembled form. Each processor 554b4319ba4SBarry Smith assembles only its local Neumann problem and the parallel matrix vector 555b4319ba4SBarry Smith product is handled "implicitly". 556b4319ba4SBarry Smith 557b4319ba4SBarry Smith Operations Provided: 5586726f965SBarry Smith + MatMult() 5592e74eeadSLisandro Dalcin . MatMultAdd() 5602e74eeadSLisandro Dalcin . MatMultTranspose() 5612e74eeadSLisandro Dalcin . MatMultTransposeAdd() 5626726f965SBarry Smith . MatZeroEntries() 5636726f965SBarry Smith . MatSetOption() 5642e74eeadSLisandro Dalcin . MatZeroRows() 5656726f965SBarry Smith . MatZeroRowsLocal() 5662e74eeadSLisandro Dalcin . MatSetValues() 5676726f965SBarry Smith . MatSetValuesLocal() 5682e74eeadSLisandro Dalcin . MatScale() 5692e74eeadSLisandro Dalcin . MatGetDiagonal() 5706726f965SBarry Smith - MatSetLocalToGlobalMapping() 571b4319ba4SBarry Smith 572b4319ba4SBarry Smith Options Database Keys: 573b4319ba4SBarry Smith . -mat_type is - sets the matrix type to "is" during a call to MatSetFromOptions() 574b4319ba4SBarry Smith 575b4319ba4SBarry Smith Notes: Options prefix for the inner matrix are given by -is_mat_xxx 576b4319ba4SBarry Smith 577b4319ba4SBarry Smith You must call MatSetLocalToGlobalMapping() before using this matrix type. 578b4319ba4SBarry Smith 579b4319ba4SBarry Smith You can do matrix preallocation on the local matrix after you obtain it with 580b4319ba4SBarry Smith MatISGetLocalMat() 581b4319ba4SBarry Smith 582b4319ba4SBarry Smith Level: advanced 583b4319ba4SBarry Smith 584b4319ba4SBarry Smith .seealso: PC, MatISGetLocalMat(), MatSetLocalToGlobalMapping() 585b4319ba4SBarry Smith 586b4319ba4SBarry Smith M*/ 587b4319ba4SBarry Smith 588b4319ba4SBarry Smith EXTERN_C_BEGIN 589b4319ba4SBarry Smith #undef __FUNCT__ 590b4319ba4SBarry Smith #define __FUNCT__ "MatCreate_IS" 5917087cfbeSBarry Smith PetscErrorCode MatCreate_IS(Mat A) 592b4319ba4SBarry Smith { 593dfbe8321SBarry Smith PetscErrorCode ierr; 594b4319ba4SBarry Smith Mat_IS *b; 595b4319ba4SBarry Smith 596b4319ba4SBarry Smith PetscFunctionBegin; 59738f2d2fdSLisandro Dalcin ierr = PetscNewLog(A,Mat_IS,&b);CHKERRQ(ierr); 598b4319ba4SBarry Smith A->data = (void*)b; 599b4319ba4SBarry Smith ierr = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 600b4319ba4SBarry Smith 601b4319ba4SBarry Smith A->ops->mult = MatMult_IS; 6022e74eeadSLisandro Dalcin A->ops->multadd = MatMultAdd_IS; 6032e74eeadSLisandro Dalcin A->ops->multtranspose = MatMultTranspose_IS; 6042e74eeadSLisandro Dalcin A->ops->multtransposeadd = MatMultTransposeAdd_IS; 605b4319ba4SBarry Smith A->ops->destroy = MatDestroy_IS; 606b4319ba4SBarry Smith A->ops->setlocaltoglobalmapping = MatSetLocalToGlobalMapping_IS; 6072e74eeadSLisandro Dalcin A->ops->setvalues = MatSetValues_IS; 608b4319ba4SBarry Smith A->ops->setvalueslocal = MatSetValuesLocal_IS; 6092e74eeadSLisandro Dalcin A->ops->zerorows = MatZeroRows_IS; 610b4319ba4SBarry Smith A->ops->zerorowslocal = MatZeroRowsLocal_IS; 611b4319ba4SBarry Smith A->ops->assemblybegin = MatAssemblyBegin_IS; 612b4319ba4SBarry Smith A->ops->assemblyend = MatAssemblyEnd_IS; 613b4319ba4SBarry Smith A->ops->view = MatView_IS; 6146726f965SBarry Smith A->ops->zeroentries = MatZeroEntries_IS; 6152e74eeadSLisandro Dalcin A->ops->scale = MatScale_IS; 6162e74eeadSLisandro Dalcin A->ops->getdiagonal = MatGetDiagonal_IS; 6176726f965SBarry Smith A->ops->setoption = MatSetOption_IS; 618*4e4c7dbeSStefano Zampini A->ops->getvecs = MatGetVecs_IS; 619b4319ba4SBarry Smith 62026283091SBarry Smith ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 62126283091SBarry Smith ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 622b4319ba4SBarry Smith 623b4319ba4SBarry Smith b->A = 0; 624b4319ba4SBarry Smith b->ctx = 0; 625b4319ba4SBarry Smith b->x = 0; 626b4319ba4SBarry Smith b->y = 0; 627b4319ba4SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatISGetLocalMat_C","MatISGetLocalMat_IS",MatISGetLocalMat_IS);CHKERRQ(ierr); 6283b03a366Sstefano_zampini ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatISSetLocalMat_C","MatISSetLocalMat_IS",MatISSetLocalMat_IS);CHKERRQ(ierr); 62917667f90SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)A,MATIS);CHKERRQ(ierr); 630b4319ba4SBarry Smith 631b4319ba4SBarry Smith PetscFunctionReturn(0); 632b4319ba4SBarry Smith } 633b4319ba4SBarry Smith EXTERN_C_END 634