1be1d678aSKris Buschelman #define PETSCMAT_DLL 2be1d678aSKris Buschelman 3b4319ba4SBarry Smith /* 4b4319ba4SBarry Smith Creates a matrix class for using the Neumann-Neumann type preconditioners. 5b4319ba4SBarry Smith This stores the matrices in globally unassembled form. Each processor 6b4319ba4SBarry Smith assembles only its local Neumann problem and the parallel matrix vector 7b4319ba4SBarry Smith product is handled "implicitly". 8b4319ba4SBarry Smith 9b4319ba4SBarry Smith We provide: 10b4319ba4SBarry Smith MatMult() 11b4319ba4SBarry Smith 12b4319ba4SBarry Smith Currently this allows for only one subdomain per processor. 13b4319ba4SBarry Smith 14b4319ba4SBarry Smith */ 15b4319ba4SBarry Smith 167c4f633dSBarry Smith #include "../src/mat/impls/is/matis.h" /*I "petscmat.h" I*/ 17b4319ba4SBarry Smith 18b4319ba4SBarry Smith #undef __FUNCT__ 19b4319ba4SBarry Smith #define __FUNCT__ "MatDestroy_IS" 20dfbe8321SBarry Smith PetscErrorCode MatDestroy_IS(Mat A) 21b4319ba4SBarry Smith { 22dfbe8321SBarry Smith PetscErrorCode ierr; 23b4319ba4SBarry Smith Mat_IS *b = (Mat_IS*)A->data; 24b4319ba4SBarry Smith 25b4319ba4SBarry Smith PetscFunctionBegin; 26b4319ba4SBarry Smith if (b->A) { 27b4319ba4SBarry Smith ierr = MatDestroy(b->A);CHKERRQ(ierr); 28b4319ba4SBarry Smith } 29b4319ba4SBarry Smith if (b->ctx) { 30b4319ba4SBarry Smith ierr = VecScatterDestroy(b->ctx);CHKERRQ(ierr); 31b4319ba4SBarry Smith } 32b4319ba4SBarry Smith if (b->x) { 33b4319ba4SBarry Smith ierr = VecDestroy(b->x);CHKERRQ(ierr); 34b4319ba4SBarry Smith } 35b4319ba4SBarry Smith if (b->y) { 36b4319ba4SBarry Smith ierr = VecDestroy(b->y);CHKERRQ(ierr); 37b4319ba4SBarry Smith } 38b4319ba4SBarry Smith if (b->mapping) { 39b4319ba4SBarry Smith ierr = ISLocalToGlobalMappingDestroy(b->mapping);CHKERRQ(ierr); 40b4319ba4SBarry Smith } 41b4319ba4SBarry Smith ierr = PetscFree(b);CHKERRQ(ierr); 42dbd8c25aSHong Zhang ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr); 43901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C","",PETSC_NULL);CHKERRQ(ierr); 44b4319ba4SBarry Smith PetscFunctionReturn(0); 45b4319ba4SBarry Smith } 46b4319ba4SBarry Smith 47b4319ba4SBarry Smith #undef __FUNCT__ 48b4319ba4SBarry Smith #define __FUNCT__ "MatMult_IS" 49dfbe8321SBarry Smith PetscErrorCode MatMult_IS(Mat A,Vec x,Vec y) 50b4319ba4SBarry Smith { 51dfbe8321SBarry Smith PetscErrorCode ierr; 52b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 53b4319ba4SBarry Smith PetscScalar zero = 0.0; 54b4319ba4SBarry Smith 55b4319ba4SBarry Smith PetscFunctionBegin; 56b4319ba4SBarry Smith /* scatter the global vector x into the local work vector */ 57ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 58ca9f406cSSatish Balay ierr = VecScatterEnd(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 59b4319ba4SBarry Smith 60b4319ba4SBarry Smith /* multiply the local matrix */ 61b4319ba4SBarry Smith ierr = MatMult(is->A,is->x,is->y);CHKERRQ(ierr); 62b4319ba4SBarry Smith 63b4319ba4SBarry Smith /* scatter product back into global memory */ 642dcb1b2aSMatthew Knepley ierr = VecSet(y,zero);CHKERRQ(ierr); 65ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 66ca9f406cSSatish Balay ierr = VecScatterEnd(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 67b4319ba4SBarry Smith 68b4319ba4SBarry Smith PetscFunctionReturn(0); 69b4319ba4SBarry Smith } 70b4319ba4SBarry Smith 71b4319ba4SBarry Smith #undef __FUNCT__ 722e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultAdd_IS" 732e74eeadSLisandro Dalcin static PetscErrorCode MatMultAdd_IS(Mat A,Vec v1,Vec v2,Vec v3) 742e74eeadSLisandro Dalcin { 752e74eeadSLisandro Dalcin PetscErrorCode ierr; 762e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 772e74eeadSLisandro Dalcin 782e74eeadSLisandro Dalcin PetscFunctionBegin; /* v3 = v2 + A * v1.*/ 792e74eeadSLisandro Dalcin /* scatter the global vector v1 into the local work vector */ 80ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,v1,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 81ca9f406cSSatish Balay ierr = VecScatterEnd (is->ctx,v1,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 82ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,v2,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 83ca9f406cSSatish Balay ierr = VecScatterEnd (is->ctx,v2,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 842e74eeadSLisandro Dalcin 852e74eeadSLisandro Dalcin /* multiply the local matrix */ 862e74eeadSLisandro Dalcin ierr = MatMultAdd(is->A,is->x,is->y,is->y);CHKERRQ(ierr); 872e74eeadSLisandro Dalcin 882e74eeadSLisandro Dalcin /* scatter result back into global vector */ 892e74eeadSLisandro Dalcin ierr = VecSet(v3,0);CHKERRQ(ierr); 90ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,is->y,v3,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 91ca9f406cSSatish Balay ierr = VecScatterEnd (is->ctx,is->y,v3,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 922e74eeadSLisandro Dalcin PetscFunctionReturn(0); 932e74eeadSLisandro Dalcin } 942e74eeadSLisandro Dalcin 952e74eeadSLisandro Dalcin #undef __FUNCT__ 962e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTranspose_IS" 972e74eeadSLisandro Dalcin PetscErrorCode MatMultTranspose_IS(Mat A,Vec x,Vec y) 982e74eeadSLisandro Dalcin { 992e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 1002e74eeadSLisandro Dalcin PetscErrorCode ierr; 1012e74eeadSLisandro Dalcin 1022e74eeadSLisandro Dalcin PetscFunctionBegin; /* y = A' * x */ 1032e74eeadSLisandro Dalcin /* scatter the global vector x into the local work vector */ 104ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 105ca9f406cSSatish Balay ierr = VecScatterEnd(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1062e74eeadSLisandro Dalcin 1072e74eeadSLisandro Dalcin /* multiply the local matrix */ 1082e74eeadSLisandro Dalcin ierr = MatMultTranspose(is->A,is->x,is->y);CHKERRQ(ierr); 1092e74eeadSLisandro Dalcin 1102e74eeadSLisandro Dalcin /* scatter product back into global vector */ 1112e74eeadSLisandro Dalcin ierr = VecSet(y,0);CHKERRQ(ierr); 112ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 113ca9f406cSSatish Balay ierr = VecScatterEnd(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1142e74eeadSLisandro Dalcin PetscFunctionReturn(0); 1152e74eeadSLisandro Dalcin } 1162e74eeadSLisandro Dalcin 1172e74eeadSLisandro Dalcin #undef __FUNCT__ 1182e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTransposeAdd_IS" 1192e74eeadSLisandro Dalcin PetscErrorCode MatMultTransposeAdd_IS(Mat A,Vec v1,Vec v2,Vec v3) 1202e74eeadSLisandro Dalcin { 1212e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 1222e74eeadSLisandro Dalcin PetscErrorCode ierr; 1232e74eeadSLisandro Dalcin 1242e74eeadSLisandro Dalcin PetscFunctionBegin; /* v3 = v2 + A' * v1.*/ 1252e74eeadSLisandro Dalcin /* scatter the global vectors v1 and v2 into the local work vectors */ 126ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,v1,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 127ca9f406cSSatish Balay ierr = VecScatterEnd (is->ctx,v1,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 128ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,v2,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 129ca9f406cSSatish Balay ierr = VecScatterEnd (is->ctx,v2,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1302e74eeadSLisandro Dalcin 1312e74eeadSLisandro Dalcin /* multiply the local matrix */ 1322e74eeadSLisandro Dalcin ierr = MatMultTransposeAdd(is->A,is->x,is->y,is->y);CHKERRQ(ierr); 1332e74eeadSLisandro Dalcin 1342e74eeadSLisandro Dalcin /* scatter result back into global vector */ 1352e74eeadSLisandro Dalcin ierr = VecSet(v3,0);CHKERRQ(ierr); 136ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,is->y,v3,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 137ca9f406cSSatish Balay ierr = VecScatterEnd (is->ctx,is->y,v3,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1382e74eeadSLisandro Dalcin PetscFunctionReturn(0); 1392e74eeadSLisandro Dalcin } 1402e74eeadSLisandro Dalcin 1412e74eeadSLisandro Dalcin #undef __FUNCT__ 142b4319ba4SBarry Smith #define __FUNCT__ "MatView_IS" 143dfbe8321SBarry Smith PetscErrorCode MatView_IS(Mat A,PetscViewer viewer) 144b4319ba4SBarry Smith { 145b4319ba4SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 146dfbe8321SBarry Smith PetscErrorCode ierr; 147b4319ba4SBarry Smith PetscViewer sviewer; 148b4319ba4SBarry Smith 149b4319ba4SBarry Smith PetscFunctionBegin; 150b4319ba4SBarry Smith ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr); 151b4319ba4SBarry Smith ierr = MatView(a->A,sviewer);CHKERRQ(ierr); 152b4319ba4SBarry Smith ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr); 153b4319ba4SBarry Smith PetscFunctionReturn(0); 154b4319ba4SBarry Smith } 155b4319ba4SBarry Smith 156b4319ba4SBarry Smith #undef __FUNCT__ 157b4319ba4SBarry Smith #define __FUNCT__ "MatSetLocalToGlobalMapping_IS" 158dfbe8321SBarry Smith PetscErrorCode MatSetLocalToGlobalMapping_IS(Mat A,ISLocalToGlobalMapping mapping) 159b4319ba4SBarry Smith { 160dfbe8321SBarry Smith PetscErrorCode ierr; 16113f74950SBarry Smith PetscInt n; 162b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 163b4319ba4SBarry Smith IS from,to; 164b4319ba4SBarry Smith Vec global; 165b4319ba4SBarry Smith 166b4319ba4SBarry Smith PetscFunctionBegin; 167e7e72b3dSBarry Smith if (is->mapping) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Mapping already set for matrix"); 1682e74eeadSLisandro Dalcin PetscCheckSameComm(A,1,mapping,2); 169b4319ba4SBarry Smith ierr = PetscObjectReference((PetscObject)mapping);CHKERRQ(ierr); 170c3122656SLisandro Dalcin if (is->mapping) { ierr = ISLocalToGlobalMappingDestroy(is->mapping);CHKERRQ(ierr); } 171c3122656SLisandro Dalcin is->mapping = mapping; 172b4319ba4SBarry Smith 173b4319ba4SBarry Smith /* Create the local matrix A */ 174b4319ba4SBarry Smith ierr = ISLocalToGlobalMappingGetSize(mapping,&n);CHKERRQ(ierr); 175f69a0ea3SMatthew Knepley ierr = MatCreate(PETSC_COMM_SELF,&is->A);CHKERRQ(ierr); 176f69a0ea3SMatthew Knepley ierr = MatSetSizes(is->A,n,n,n,n);CHKERRQ(ierr); 1772e74eeadSLisandro Dalcin ierr = MatSetOptionsPrefix(is->A,"is");CHKERRQ(ierr); 178b4319ba4SBarry Smith ierr = MatSetFromOptions(is->A);CHKERRQ(ierr); 179b4319ba4SBarry Smith 180b4319ba4SBarry Smith /* Create the local work vectors */ 181b4319ba4SBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,n,&is->x);CHKERRQ(ierr); 182b4319ba4SBarry Smith ierr = VecDuplicate(is->x,&is->y);CHKERRQ(ierr); 183b4319ba4SBarry Smith 184b4319ba4SBarry Smith /* setup the global to local scatter */ 185b4319ba4SBarry Smith ierr = ISCreateStride(PETSC_COMM_SELF,n,0,1,&to);CHKERRQ(ierr); 186b4319ba4SBarry Smith ierr = ISLocalToGlobalMappingApplyIS(mapping,to,&from);CHKERRQ(ierr); 187d0f46423SBarry Smith ierr = VecCreateMPIWithArray(((PetscObject)A)->comm,A->cmap->n,A->cmap->N,PETSC_NULL,&global);CHKERRQ(ierr); 188b4319ba4SBarry Smith ierr = VecScatterCreate(global,from,is->x,to,&is->ctx);CHKERRQ(ierr); 189b4319ba4SBarry Smith ierr = VecDestroy(global);CHKERRQ(ierr); 190b4319ba4SBarry Smith ierr = ISDestroy(to);CHKERRQ(ierr); 191b4319ba4SBarry Smith ierr = ISDestroy(from);CHKERRQ(ierr); 192b4319ba4SBarry Smith PetscFunctionReturn(0); 193b4319ba4SBarry Smith } 194b4319ba4SBarry Smith 1952e74eeadSLisandro Dalcin #define ISG2LMapApply(mapping,n,in,out) 0;\ 1962e74eeadSLisandro Dalcin if (!(mapping)->globals) {\ 1972e74eeadSLisandro Dalcin PetscErrorCode _ierr = ISGlobalToLocalMappingApply((mapping),IS_GTOLM_MASK,0,0,0,0);CHKERRQ(_ierr);\ 1982e74eeadSLisandro Dalcin }\ 1992e74eeadSLisandro Dalcin {\ 2002e74eeadSLisandro Dalcin PetscInt _i,*_globals = (mapping)->globals,_start = (mapping)->globalstart,_end = (mapping)->globalend;\ 2012e74eeadSLisandro Dalcin for (_i=0; _i<n; _i++) {\ 2022e74eeadSLisandro Dalcin if (in[_i] < 0) out[_i] = in[_i];\ 2032e74eeadSLisandro Dalcin else if (in[_i] < _start) out[_i] = -1;\ 2042e74eeadSLisandro Dalcin else if (in[_i] > _end) out[_i] = -1;\ 2052e74eeadSLisandro Dalcin else out[_i] = _globals[in[_i] - _start];\ 2062e74eeadSLisandro Dalcin }\ 2072e74eeadSLisandro Dalcin } 2082e74eeadSLisandro Dalcin 2092e74eeadSLisandro Dalcin #undef __FUNCT__ 2102e74eeadSLisandro Dalcin #define __FUNCT__ "MatSetValues_IS" 2112e74eeadSLisandro Dalcin PetscErrorCode MatSetValues_IS(Mat mat, PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols, const PetscScalar *values, InsertMode addv) 2122e74eeadSLisandro Dalcin { 2132e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)mat->data; 2142e74eeadSLisandro Dalcin PetscInt rows_l[2048],cols_l[2048]; 2152e74eeadSLisandro Dalcin PetscErrorCode ierr; 2162e74eeadSLisandro Dalcin 2172e74eeadSLisandro Dalcin PetscFunctionBegin; 2182e74eeadSLisandro Dalcin #if defined(PETSC_USE_DEBUG) 2192e74eeadSLisandro Dalcin if (m > 2048 || n > 2048) { 220e32f2f54SBarry Smith SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Number of row/column indices must be <= 2048: they are %D %D",m,n); 2212e74eeadSLisandro Dalcin } 2222e74eeadSLisandro Dalcin #endif 2232e74eeadSLisandro Dalcin ierr = ISG2LMapApply(is->mapping,m,rows,rows_l);CHKERRQ(ierr); 2242e74eeadSLisandro Dalcin ierr = ISG2LMapApply(is->mapping,n,cols,cols_l);CHKERRQ(ierr); 2252e74eeadSLisandro Dalcin ierr = MatSetValues(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr); 2262e74eeadSLisandro Dalcin PetscFunctionReturn(0); 2272e74eeadSLisandro Dalcin } 2282e74eeadSLisandro Dalcin 2292e74eeadSLisandro Dalcin #undef ISG2LMapSetUp 2302e74eeadSLisandro Dalcin #undef ISG2LMapApply 231b4319ba4SBarry Smith 232b4319ba4SBarry Smith #undef __FUNCT__ 233b4319ba4SBarry Smith #define __FUNCT__ "MatSetValuesLocal_IS" 23413f74950SBarry Smith PetscErrorCode MatSetValuesLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv) 235b4319ba4SBarry Smith { 236dfbe8321SBarry Smith PetscErrorCode ierr; 237b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 238b4319ba4SBarry Smith 239b4319ba4SBarry Smith PetscFunctionBegin; 240b4319ba4SBarry Smith ierr = MatSetValues(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr); 241b4319ba4SBarry Smith PetscFunctionReturn(0); 242b4319ba4SBarry Smith } 243b4319ba4SBarry Smith 244b4319ba4SBarry Smith #undef __FUNCT__ 2452e74eeadSLisandro Dalcin #define __FUNCT__ "MatZeroRows_IS" 2462b40b63fSBarry Smith PetscErrorCode MatZeroRows_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b) 2472e74eeadSLisandro Dalcin { 2482e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 2492e74eeadSLisandro Dalcin PetscInt n_l=0, *rows_l = PETSC_NULL; 2502e74eeadSLisandro Dalcin PetscErrorCode ierr; 2512e74eeadSLisandro Dalcin 2522e74eeadSLisandro Dalcin PetscFunctionBegin; 253*97b48c8fSBarry Smith if (x && b) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"No support"); 2542e74eeadSLisandro Dalcin if (n) { 2552e74eeadSLisandro Dalcin ierr = PetscMalloc(n*sizeof(PetscInt),&rows_l);CHKERRQ(ierr); 2562e74eeadSLisandro Dalcin ierr = ISGlobalToLocalMappingApply(is->mapping,IS_GTOLM_DROP,n,rows,&n_l,rows_l);CHKERRQ(ierr); 2572e74eeadSLisandro Dalcin } 2582b40b63fSBarry Smith ierr = MatZeroRowsLocal(A,n_l,rows_l,diag,x,b);CHKERRQ(ierr); 2592e74eeadSLisandro Dalcin if (rows_l) { ierr = PetscFree(rows_l);CHKERRQ(ierr); } 2602e74eeadSLisandro Dalcin PetscFunctionReturn(0); 2612e74eeadSLisandro Dalcin } 2622e74eeadSLisandro Dalcin 2632e74eeadSLisandro Dalcin #undef __FUNCT__ 264b4319ba4SBarry Smith #define __FUNCT__ "MatZeroRowsLocal_IS" 2652b40b63fSBarry Smith PetscErrorCode MatZeroRowsLocal_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b) 266b4319ba4SBarry Smith { 267b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 268dfbe8321SBarry Smith PetscErrorCode ierr; 269f4df32b1SMatthew Knepley PetscInt i; 270b4319ba4SBarry Smith PetscScalar *array; 271b4319ba4SBarry Smith 272b4319ba4SBarry Smith PetscFunctionBegin; 273b4319ba4SBarry Smith { 274b4319ba4SBarry Smith /* 275b4319ba4SBarry Smith Set up is->x as a "counting vector". This is in order to MatMult_IS 276b4319ba4SBarry Smith work properly in the interface nodes. 277b4319ba4SBarry Smith */ 278b4319ba4SBarry Smith Vec counter; 279b4319ba4SBarry Smith PetscScalar one=1.0, zero=0.0; 280d0f46423SBarry Smith ierr = VecCreateMPI(((PetscObject)A)->comm,A->cmap->n,A->cmap->N,&counter);CHKERRQ(ierr); 2812dcb1b2aSMatthew Knepley ierr = VecSet(counter,zero);CHKERRQ(ierr); 2822dcb1b2aSMatthew Knepley ierr = VecSet(is->x,one);CHKERRQ(ierr); 283ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,is->x,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 284ca9f406cSSatish Balay ierr = VecScatterEnd (is->ctx,is->x,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 285ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,counter,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 286ca9f406cSSatish Balay ierr = VecScatterEnd (is->ctx,counter,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 287b4319ba4SBarry Smith ierr = VecDestroy(counter);CHKERRQ(ierr); 288b4319ba4SBarry Smith } 289958c9bccSBarry Smith if (!n) { 290b4319ba4SBarry Smith is->pure_neumann = PETSC_TRUE; 291b4319ba4SBarry Smith } else { 292b4319ba4SBarry Smith is->pure_neumann = PETSC_FALSE; 293b4319ba4SBarry Smith ierr = VecGetArray(is->x,&array);CHKERRQ(ierr); 2942b40b63fSBarry Smith ierr = MatZeroRows(is->A,n,rows,diag,0,0);CHKERRQ(ierr); 295b4319ba4SBarry Smith for (i=0; i<n; i++) { 296f4df32b1SMatthew Knepley ierr = MatSetValue(is->A,rows[i],rows[i],diag/(array[rows[i]]),INSERT_VALUES);CHKERRQ(ierr); 297b4319ba4SBarry Smith } 298b4319ba4SBarry Smith ierr = MatAssemblyBegin(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 299b4319ba4SBarry Smith ierr = MatAssemblyEnd (is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 300b4319ba4SBarry Smith ierr = VecRestoreArray(is->x,&array);CHKERRQ(ierr); 301b4319ba4SBarry Smith } 302b4319ba4SBarry Smith 303b4319ba4SBarry Smith PetscFunctionReturn(0); 304b4319ba4SBarry Smith } 305b4319ba4SBarry Smith 306b4319ba4SBarry Smith #undef __FUNCT__ 307b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyBegin_IS" 308dfbe8321SBarry Smith PetscErrorCode MatAssemblyBegin_IS(Mat A,MatAssemblyType type) 309b4319ba4SBarry Smith { 310b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 311dfbe8321SBarry Smith PetscErrorCode ierr; 312b4319ba4SBarry Smith 313b4319ba4SBarry Smith PetscFunctionBegin; 314b4319ba4SBarry Smith ierr = MatAssemblyBegin(is->A,type);CHKERRQ(ierr); 315b4319ba4SBarry Smith PetscFunctionReturn(0); 316b4319ba4SBarry Smith } 317b4319ba4SBarry Smith 318b4319ba4SBarry Smith #undef __FUNCT__ 319b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyEnd_IS" 320dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_IS(Mat A,MatAssemblyType type) 321b4319ba4SBarry Smith { 322b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 323dfbe8321SBarry Smith PetscErrorCode ierr; 324b4319ba4SBarry Smith 325b4319ba4SBarry Smith PetscFunctionBegin; 326b4319ba4SBarry Smith ierr = MatAssemblyEnd(is->A,type);CHKERRQ(ierr); 327b4319ba4SBarry Smith PetscFunctionReturn(0); 328b4319ba4SBarry Smith } 329b4319ba4SBarry Smith 330b4319ba4SBarry Smith EXTERN_C_BEGIN 331b4319ba4SBarry Smith #undef __FUNCT__ 332b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat_IS" 333be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatISGetLocalMat_IS(Mat mat,Mat *local) 334b4319ba4SBarry Smith { 335b4319ba4SBarry Smith Mat_IS *is = (Mat_IS *)mat->data; 336b4319ba4SBarry Smith 337b4319ba4SBarry Smith PetscFunctionBegin; 338b4319ba4SBarry Smith *local = is->A; 339b4319ba4SBarry Smith PetscFunctionReturn(0); 340b4319ba4SBarry Smith } 341b4319ba4SBarry Smith EXTERN_C_END 342b4319ba4SBarry Smith 343b4319ba4SBarry Smith #undef __FUNCT__ 344b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat" 345b4319ba4SBarry Smith /*@ 346b4319ba4SBarry Smith MatISGetLocalMat - Gets the local matrix stored inside a MATIS matrix. 347b4319ba4SBarry Smith 348b4319ba4SBarry Smith Input Parameter: 349b4319ba4SBarry Smith . mat - the matrix 350b4319ba4SBarry Smith 351b4319ba4SBarry Smith Output Parameter: 352b4319ba4SBarry Smith . local - the local matrix usually MATSEQAIJ 353b4319ba4SBarry Smith 354b4319ba4SBarry Smith Level: advanced 355b4319ba4SBarry Smith 356b4319ba4SBarry Smith Notes: 357b4319ba4SBarry Smith This can be called if you have precomputed the nonzero structure of the 358b4319ba4SBarry Smith matrix and want to provide it to the inner matrix object to improve the performance 359b4319ba4SBarry Smith of the MatSetValues() operation. 360b4319ba4SBarry Smith 361b4319ba4SBarry Smith .seealso: MATIS 362b4319ba4SBarry Smith @*/ 363be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatISGetLocalMat(Mat mat,Mat *local) 364b4319ba4SBarry Smith { 3654ac538c5SBarry Smith PetscErrorCode ierr; 366b4319ba4SBarry Smith 367b4319ba4SBarry Smith PetscFunctionBegin; 3680700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 369b4319ba4SBarry Smith PetscValidPointer(local,2); 3704ac538c5SBarry Smith ierr = PetscUseMethod(mat,"MatISGetLocalMat_C",(Mat,Mat *),(mat,local));CHKERRQ(ierr); 371b4319ba4SBarry Smith PetscFunctionReturn(0); 372b4319ba4SBarry Smith } 373b4319ba4SBarry Smith 3746726f965SBarry Smith #undef __FUNCT__ 3756726f965SBarry Smith #define __FUNCT__ "MatZeroEntries_IS" 3766726f965SBarry Smith PetscErrorCode MatZeroEntries_IS(Mat A) 3776726f965SBarry Smith { 3786726f965SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 3796726f965SBarry Smith PetscErrorCode ierr; 3806726f965SBarry Smith 3816726f965SBarry Smith PetscFunctionBegin; 3826726f965SBarry Smith ierr = MatZeroEntries(a->A);CHKERRQ(ierr); 3836726f965SBarry Smith PetscFunctionReturn(0); 3846726f965SBarry Smith } 3856726f965SBarry Smith 3866726f965SBarry Smith #undef __FUNCT__ 3872e74eeadSLisandro Dalcin #define __FUNCT__ "MatScale_IS" 3882e74eeadSLisandro Dalcin PetscErrorCode MatScale_IS(Mat A,PetscScalar a) 3892e74eeadSLisandro Dalcin { 3902e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 3912e74eeadSLisandro Dalcin PetscErrorCode ierr; 3922e74eeadSLisandro Dalcin 3932e74eeadSLisandro Dalcin PetscFunctionBegin; 3942e74eeadSLisandro Dalcin ierr = MatScale(is->A,a);CHKERRQ(ierr); 3952e74eeadSLisandro Dalcin PetscFunctionReturn(0); 3962e74eeadSLisandro Dalcin } 3972e74eeadSLisandro Dalcin 3982e74eeadSLisandro Dalcin #undef __FUNCT__ 3992e74eeadSLisandro Dalcin #define __FUNCT__ "MatGetDiagonal_IS" 4002e74eeadSLisandro Dalcin PetscErrorCode MatGetDiagonal_IS(Mat A, Vec v) 4012e74eeadSLisandro Dalcin { 4022e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 4032e74eeadSLisandro Dalcin PetscErrorCode ierr; 4042e74eeadSLisandro Dalcin 4052e74eeadSLisandro Dalcin PetscFunctionBegin; 4062e74eeadSLisandro Dalcin /* get diagonal of the local matrix */ 4072e74eeadSLisandro Dalcin ierr = MatGetDiagonal(is->A,is->x);CHKERRQ(ierr); 4082e74eeadSLisandro Dalcin 4092e74eeadSLisandro Dalcin /* scatter diagonal back into global vector */ 4102e74eeadSLisandro Dalcin ierr = VecSet(v,0);CHKERRQ(ierr); 411ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,is->x,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 412ca9f406cSSatish Balay ierr = VecScatterEnd (is->ctx,is->x,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 4132e74eeadSLisandro Dalcin PetscFunctionReturn(0); 4142e74eeadSLisandro Dalcin } 4152e74eeadSLisandro Dalcin 4162e74eeadSLisandro Dalcin #undef __FUNCT__ 4176726f965SBarry Smith #define __FUNCT__ "MatSetOption_IS" 418ace3abfcSBarry Smith PetscErrorCode MatSetOption_IS(Mat A,MatOption op,PetscBool flg) 4196726f965SBarry Smith { 4206726f965SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 4216726f965SBarry Smith PetscErrorCode ierr; 4226726f965SBarry Smith 4236726f965SBarry Smith PetscFunctionBegin; 4244e0d8c25SBarry Smith ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr); 4256726f965SBarry Smith PetscFunctionReturn(0); 4266726f965SBarry Smith } 4276726f965SBarry Smith 428284134d9SBarry Smith #undef __FUNCT__ 429284134d9SBarry Smith #define __FUNCT__ "MatCreateIS" 430284134d9SBarry Smith /*@ 431284134d9SBarry Smith MatCreateIS - Creates a "process" unassmembled matrix, it is assembled on each 432284134d9SBarry Smith process but not across processes. 433284134d9SBarry Smith 434284134d9SBarry Smith Input Parameters: 435284134d9SBarry Smith + comm - MPI communicator that will share the matrix 436284134d9SBarry Smith . m,n,M,N - local and/or global sizes of the the left and right vector used in matrix vector products 437284134d9SBarry Smith - map - mapping that defines the global number for each local number 438284134d9SBarry Smith 439284134d9SBarry Smith Output Parameter: 440284134d9SBarry Smith . A - the resulting matrix 441284134d9SBarry Smith 4428e6c10adSSatish Balay Level: advanced 4438e6c10adSSatish Balay 444284134d9SBarry Smith Notes: See MATIS for more details 4458cda6cd7SBarry Smith m and n are NOT related to the size of the map, they are the size of the part of the vector owned 4468cda6cd7SBarry Smith by that process. m + nghosts (or n + nghosts) is the length of map since map maps all local points 4478cda6cd7SBarry Smith plus the ghost points to global indices. 448284134d9SBarry Smith 449284134d9SBarry Smith .seealso: MATIS, MatSetLocalToGlobalMapping() 450284134d9SBarry Smith @*/ 451284134d9SBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatCreateIS(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,ISLocalToGlobalMapping map,Mat *A) 452284134d9SBarry Smith { 453284134d9SBarry Smith PetscErrorCode ierr; 454284134d9SBarry Smith 455284134d9SBarry Smith PetscFunctionBegin; 456284134d9SBarry Smith ierr = MatCreate(comm,A);CHKERRQ(ierr); 457284134d9SBarry Smith ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr); 458284134d9SBarry Smith ierr = MatSetType(*A,MATIS);CHKERRQ(ierr); 459284134d9SBarry Smith ierr = MatSetLocalToGlobalMapping(*A,map);CHKERRQ(ierr); 460284134d9SBarry Smith PetscFunctionReturn(0); 461284134d9SBarry Smith } 462284134d9SBarry Smith 463b4319ba4SBarry Smith /*MC 464b4319ba4SBarry Smith MATIS - MATIS = "is" - A matrix type to be used for using the Neumann-Neumann type preconditioners. 465b4319ba4SBarry Smith This stores the matrices in globally unassembled form. Each processor 466b4319ba4SBarry Smith assembles only its local Neumann problem and the parallel matrix vector 467b4319ba4SBarry Smith product is handled "implicitly". 468b4319ba4SBarry Smith 469b4319ba4SBarry Smith Operations Provided: 4706726f965SBarry Smith + MatMult() 4712e74eeadSLisandro Dalcin . MatMultAdd() 4722e74eeadSLisandro Dalcin . MatMultTranspose() 4732e74eeadSLisandro Dalcin . MatMultTransposeAdd() 4746726f965SBarry Smith . MatZeroEntries() 4756726f965SBarry Smith . MatSetOption() 4762e74eeadSLisandro Dalcin . MatZeroRows() 4776726f965SBarry Smith . MatZeroRowsLocal() 4782e74eeadSLisandro Dalcin . MatSetValues() 4796726f965SBarry Smith . MatSetValuesLocal() 4802e74eeadSLisandro Dalcin . MatScale() 4812e74eeadSLisandro Dalcin . MatGetDiagonal() 4826726f965SBarry Smith - MatSetLocalToGlobalMapping() 483b4319ba4SBarry Smith 484b4319ba4SBarry Smith Options Database Keys: 485b4319ba4SBarry Smith . -mat_type is - sets the matrix type to "is" during a call to MatSetFromOptions() 486b4319ba4SBarry Smith 487b4319ba4SBarry Smith Notes: Options prefix for the inner matrix are given by -is_mat_xxx 488b4319ba4SBarry Smith 489b4319ba4SBarry Smith You must call MatSetLocalToGlobalMapping() before using this matrix type. 490b4319ba4SBarry Smith 491b4319ba4SBarry Smith You can do matrix preallocation on the local matrix after you obtain it with 492b4319ba4SBarry Smith MatISGetLocalMat() 493b4319ba4SBarry Smith 494b4319ba4SBarry Smith Level: advanced 495b4319ba4SBarry Smith 496b4319ba4SBarry Smith .seealso: PC, MatISGetLocalMat(), MatSetLocalToGlobalMapping() 497b4319ba4SBarry Smith 498b4319ba4SBarry Smith M*/ 499b4319ba4SBarry Smith 500b4319ba4SBarry Smith EXTERN_C_BEGIN 501b4319ba4SBarry Smith #undef __FUNCT__ 502b4319ba4SBarry Smith #define __FUNCT__ "MatCreate_IS" 503be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_IS(Mat A) 504b4319ba4SBarry Smith { 505dfbe8321SBarry Smith PetscErrorCode ierr; 506b4319ba4SBarry Smith Mat_IS *b; 507b4319ba4SBarry Smith 508b4319ba4SBarry Smith PetscFunctionBegin; 50938f2d2fdSLisandro Dalcin ierr = PetscNewLog(A,Mat_IS,&b);CHKERRQ(ierr); 510b4319ba4SBarry Smith A->data = (void*)b; 511b4319ba4SBarry Smith ierr = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 512b4319ba4SBarry Smith A->mapping = 0; 513b4319ba4SBarry Smith 514b4319ba4SBarry Smith A->ops->mult = MatMult_IS; 5152e74eeadSLisandro Dalcin A->ops->multadd = MatMultAdd_IS; 5162e74eeadSLisandro Dalcin A->ops->multtranspose = MatMultTranspose_IS; 5172e74eeadSLisandro Dalcin A->ops->multtransposeadd = MatMultTransposeAdd_IS; 518b4319ba4SBarry Smith A->ops->destroy = MatDestroy_IS; 519b4319ba4SBarry Smith A->ops->setlocaltoglobalmapping = MatSetLocalToGlobalMapping_IS; 5202e74eeadSLisandro Dalcin A->ops->setvalues = MatSetValues_IS; 521b4319ba4SBarry Smith A->ops->setvalueslocal = MatSetValuesLocal_IS; 5222e74eeadSLisandro Dalcin A->ops->zerorows = MatZeroRows_IS; 523b4319ba4SBarry Smith A->ops->zerorowslocal = MatZeroRowsLocal_IS; 524b4319ba4SBarry Smith A->ops->assemblybegin = MatAssemblyBegin_IS; 525b4319ba4SBarry Smith A->ops->assemblyend = MatAssemblyEnd_IS; 526b4319ba4SBarry Smith A->ops->view = MatView_IS; 5276726f965SBarry Smith A->ops->zeroentries = MatZeroEntries_IS; 5282e74eeadSLisandro Dalcin A->ops->scale = MatScale_IS; 5292e74eeadSLisandro Dalcin A->ops->getdiagonal = MatGetDiagonal_IS; 5306726f965SBarry Smith A->ops->setoption = MatSetOption_IS; 531b4319ba4SBarry Smith 53226283091SBarry Smith ierr = PetscLayoutSetBlockSize(A->rmap,1);CHKERRQ(ierr); 53326283091SBarry Smith ierr = PetscLayoutSetBlockSize(A->cmap,1);CHKERRQ(ierr); 53426283091SBarry Smith ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 53526283091SBarry Smith ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 536b4319ba4SBarry Smith 537b4319ba4SBarry Smith b->A = 0; 538b4319ba4SBarry Smith b->ctx = 0; 539b4319ba4SBarry Smith b->x = 0; 540b4319ba4SBarry Smith b->y = 0; 541b4319ba4SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatISGetLocalMat_C","MatISGetLocalMat_IS",MatISGetLocalMat_IS);CHKERRQ(ierr); 54217667f90SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)A,MATIS);CHKERRQ(ierr); 543b4319ba4SBarry Smith 544b4319ba4SBarry Smith PetscFunctionReturn(0); 545b4319ba4SBarry Smith } 546b4319ba4SBarry Smith EXTERN_C_END 547