1 /* 2 Creates a matrix class for using the Neumann-Neumann type preconditioners. 3 This stores the matrices in globally unassembled form. Each processor 4 assembles only its local Neumann problem and the parallel matrix vector 5 product is handled "implicitly". 6 7 We provide: 8 MatMult() 9 10 Currently this allows for only one subdomain per processor. 11 12 */ 13 14 #include "src/mat/impls/is/matis.h" /*I "petscmat.h" I*/ 15 16 #undef __FUNCT__ 17 #define __FUNCT__ "MatDestroy_IS" 18 PetscErrorCode MatDestroy_IS(Mat A) 19 { 20 PetscErrorCode ierr; 21 Mat_IS *b = (Mat_IS*)A->data; 22 23 PetscFunctionBegin; 24 if (b->A) { 25 ierr = MatDestroy(b->A);CHKERRQ(ierr); 26 } 27 if (b->ctx) { 28 ierr = VecScatterDestroy(b->ctx);CHKERRQ(ierr); 29 } 30 if (b->x) { 31 ierr = VecDestroy(b->x);CHKERRQ(ierr); 32 } 33 if (b->y) { 34 ierr = VecDestroy(b->y);CHKERRQ(ierr); 35 } 36 if (b->mapping) { 37 ierr = ISLocalToGlobalMappingDestroy(b->mapping);CHKERRQ(ierr); 38 } 39 ierr = PetscFree(b);CHKERRQ(ierr); 40 ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C","",PETSC_NULL);CHKERRQ(ierr); 41 PetscFunctionReturn(0); 42 } 43 44 #undef __FUNCT__ 45 #define __FUNCT__ "MatMult_IS" 46 PetscErrorCode MatMult_IS(Mat A,Vec x,Vec y) 47 { 48 PetscErrorCode ierr; 49 Mat_IS *is = (Mat_IS*)A->data; 50 PetscScalar zero = 0.0; 51 52 PetscFunctionBegin; 53 /* scatter the global vector x into the local work vector */ 54 ierr = VecScatterBegin(x,is->x,INSERT_VALUES,SCATTER_FORWARD,is->ctx);CHKERRQ(ierr); 55 ierr = VecScatterEnd(x,is->x,INSERT_VALUES,SCATTER_FORWARD,is->ctx);CHKERRQ(ierr); 56 57 /* multiply the local matrix */ 58 ierr = MatMult(is->A,is->x,is->y);CHKERRQ(ierr); 59 60 /* scatter product back into global memory */ 61 ierr = VecSet(&zero,y);CHKERRQ(ierr); 62 ierr = VecScatterBegin(is->y,y,ADD_VALUES,SCATTER_REVERSE,is->ctx);CHKERRQ(ierr); 63 ierr = VecScatterEnd(is->y,y,ADD_VALUES,SCATTER_REVERSE,is->ctx);CHKERRQ(ierr); 64 65 PetscFunctionReturn(0); 66 } 67 68 #undef __FUNCT__ 69 #define __FUNCT__ "MatView_IS" 70 PetscErrorCode MatView_IS(Mat A,PetscViewer viewer) 71 { 72 Mat_IS *a = (Mat_IS*)A->data; 73 PetscErrorCode ierr; 74 PetscViewer sviewer; 75 76 PetscFunctionBegin; 77 ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr); 78 ierr = MatView(a->A,sviewer);CHKERRQ(ierr); 79 ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr); 80 PetscFunctionReturn(0); 81 } 82 83 #undef __FUNCT__ 84 #define __FUNCT__ "MatSetLocalToGlobalMapping_IS" 85 PetscErrorCode MatSetLocalToGlobalMapping_IS(Mat A,ISLocalToGlobalMapping mapping) 86 { 87 PetscErrorCode ierr; 88 int n; 89 Mat_IS *is = (Mat_IS*)A->data; 90 IS from,to; 91 Vec global; 92 93 PetscFunctionBegin; 94 is->mapping = mapping; 95 ierr = PetscObjectReference((PetscObject)mapping);CHKERRQ(ierr); 96 97 /* Create the local matrix A */ 98 ierr = ISLocalToGlobalMappingGetSize(mapping,&n);CHKERRQ(ierr); 99 ierr = MatCreate(PETSC_COMM_SELF,n,n,n,n,&is->A);CHKERRQ(ierr); 100 ierr = PetscObjectSetOptionsPrefix((PetscObject)is->A,"is");CHKERRQ(ierr); 101 ierr = MatSetFromOptions(is->A);CHKERRQ(ierr); 102 103 /* Create the local work vectors */ 104 ierr = VecCreateSeq(PETSC_COMM_SELF,n,&is->x);CHKERRQ(ierr); 105 ierr = VecDuplicate(is->x,&is->y);CHKERRQ(ierr); 106 107 /* setup the global to local scatter */ 108 ierr = ISCreateStride(PETSC_COMM_SELF,n,0,1,&to);CHKERRQ(ierr); 109 ierr = ISLocalToGlobalMappingApplyIS(mapping,to,&from);CHKERRQ(ierr); 110 ierr = VecCreateMPI(A->comm,A->n,A->N,&global);CHKERRQ(ierr); 111 ierr = VecScatterCreate(global,from,is->x,to,&is->ctx);CHKERRQ(ierr); 112 ierr = VecDestroy(global);CHKERRQ(ierr); 113 ierr = ISDestroy(to);CHKERRQ(ierr); 114 ierr = ISDestroy(from);CHKERRQ(ierr); 115 PetscFunctionReturn(0); 116 } 117 118 119 #undef __FUNCT__ 120 #define __FUNCT__ "MatSetValuesLocal_IS" 121 PetscErrorCode MatSetValuesLocal_IS(Mat A,int m,const int *rows, int n,const int *cols,const PetscScalar *values,InsertMode addv) 122 { 123 PetscErrorCode ierr; 124 Mat_IS *is = (Mat_IS*)A->data; 125 126 PetscFunctionBegin; 127 ierr = MatSetValues(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr); 128 PetscFunctionReturn(0); 129 } 130 131 #undef __FUNCT__ 132 #define __FUNCT__ "MatZeroRowsLocal_IS" 133 PetscErrorCode MatZeroRowsLocal_IS(Mat A,IS isrows,const PetscScalar *diag) 134 { 135 Mat_IS *is = (Mat_IS*)A->data; 136 PetscErrorCode ierr; 137 int i,n,*rows; 138 PetscScalar *array; 139 140 PetscFunctionBegin; 141 142 { 143 /* 144 Set up is->x as a "counting vector". This is in order to MatMult_IS 145 work properly in the interface nodes. 146 */ 147 Vec counter; 148 PetscScalar one=1.0, zero=0.0; 149 ierr = VecCreateMPI(A->comm,A->n,A->N,&counter);CHKERRQ(ierr); 150 ierr = VecSet(&zero,counter);CHKERRQ(ierr); 151 ierr = VecSet(&one,is->x);CHKERRQ(ierr); 152 ierr = VecScatterBegin(is->x,counter,ADD_VALUES,SCATTER_REVERSE,is->ctx);CHKERRQ(ierr); 153 ierr = VecScatterEnd (is->x,counter,ADD_VALUES,SCATTER_REVERSE,is->ctx);CHKERRQ(ierr); 154 ierr = VecScatterBegin(counter,is->x,INSERT_VALUES,SCATTER_FORWARD,is->ctx);CHKERRQ(ierr); 155 ierr = VecScatterEnd (counter,is->x,INSERT_VALUES,SCATTER_FORWARD,is->ctx);CHKERRQ(ierr); 156 ierr = VecDestroy(counter);CHKERRQ(ierr); 157 } 158 ierr = ISGetLocalSize(isrows,&n);CHKERRQ(ierr); 159 if (!n) { 160 is->pure_neumann = PETSC_TRUE; 161 } else { 162 is->pure_neumann = PETSC_FALSE; 163 ierr = ISGetIndices(isrows,&rows);CHKERRQ(ierr); 164 ierr = VecGetArray(is->x,&array);CHKERRQ(ierr); 165 ierr = MatZeroRows(is->A,isrows,diag);CHKERRQ(ierr); 166 for (i=0; i<n; i++) { 167 ierr = MatSetValue(is->A,rows[i],rows[i],(*diag)/(array[rows[i]]),INSERT_VALUES);CHKERRQ(ierr); 168 } 169 ierr = MatAssemblyBegin(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 170 ierr = MatAssemblyEnd (is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 171 ierr = VecRestoreArray(is->x,&array);CHKERRQ(ierr); 172 ierr = ISRestoreIndices(isrows,&rows);CHKERRQ(ierr); 173 } 174 175 PetscFunctionReturn(0); 176 } 177 178 #undef __FUNCT__ 179 #define __FUNCT__ "MatAssemblyBegin_IS" 180 PetscErrorCode MatAssemblyBegin_IS(Mat A,MatAssemblyType type) 181 { 182 Mat_IS *is = (Mat_IS*)A->data; 183 PetscErrorCode ierr; 184 185 PetscFunctionBegin; 186 ierr = MatAssemblyBegin(is->A,type);CHKERRQ(ierr); 187 PetscFunctionReturn(0); 188 } 189 190 #undef __FUNCT__ 191 #define __FUNCT__ "MatAssemblyEnd_IS" 192 PetscErrorCode MatAssemblyEnd_IS(Mat A,MatAssemblyType type) 193 { 194 Mat_IS *is = (Mat_IS*)A->data; 195 PetscErrorCode ierr; 196 197 PetscFunctionBegin; 198 ierr = MatAssemblyEnd(is->A,type);CHKERRQ(ierr); 199 PetscFunctionReturn(0); 200 } 201 202 EXTERN_C_BEGIN 203 #undef __FUNCT__ 204 #define __FUNCT__ "MatISGetLocalMat_IS" 205 PetscErrorCode MatISGetLocalMat_IS(Mat mat,Mat *local) 206 { 207 Mat_IS *is = (Mat_IS *)mat->data; 208 209 PetscFunctionBegin; 210 *local = is->A; 211 PetscFunctionReturn(0); 212 } 213 EXTERN_C_END 214 215 #undef __FUNCT__ 216 #define __FUNCT__ "MatISGetLocalMat" 217 /*@ 218 MatISGetLocalMat - Gets the local matrix stored inside a MATIS matrix. 219 220 Input Parameter: 221 . mat - the matrix 222 223 Output Parameter: 224 . local - the local matrix usually MATSEQAIJ 225 226 Level: advanced 227 228 Notes: 229 This can be called if you have precomputed the nonzero structure of the 230 matrix and want to provide it to the inner matrix object to improve the performance 231 of the MatSetValues() operation. 232 233 .seealso: MATIS 234 @*/ 235 PetscErrorCode MatISGetLocalMat(Mat mat,Mat *local) 236 { 237 PetscErrorCode ierr,(*f)(Mat,Mat *); 238 239 PetscFunctionBegin; 240 PetscValidHeaderSpecific(mat,MAT_COOKIE,1); 241 PetscValidPointer(local,2); 242 ierr = PetscObjectQueryFunction((PetscObject)mat,"MatISGetLocalMat_C",(void (**)(void))&f);CHKERRQ(ierr); 243 if (f) { 244 ierr = (*f)(mat,local);CHKERRQ(ierr); 245 } else { 246 local = 0; 247 } 248 PetscFunctionReturn(0); 249 } 250 251 /*MC 252 MATIS - MATIS = "is" - A matrix type to be used for using the Neumann-Neumann type preconditioners. 253 This stores the matrices in globally unassembled form. Each processor 254 assembles only its local Neumann problem and the parallel matrix vector 255 product is handled "implicitly". 256 257 Operations Provided: 258 . MatMult 259 260 Options Database Keys: 261 . -mat_type is - sets the matrix type to "is" during a call to MatSetFromOptions() 262 263 Notes: Options prefix for the inner matrix are given by -is_mat_xxx 264 265 You must call MatSetLocalToGlobalMapping() before using this matrix type. 266 267 You can do matrix preallocation on the local matrix after you obtain it with 268 MatISGetLocalMat() 269 270 Level: advanced 271 272 .seealso: PC, MatISGetLocalMat(), MatSetLocalToGlobalMapping() 273 274 M*/ 275 276 EXTERN_C_BEGIN 277 #undef __FUNCT__ 278 #define __FUNCT__ "MatCreate_IS" 279 PetscErrorCode MatCreate_IS(Mat A) 280 { 281 PetscErrorCode ierr; 282 Mat_IS *b; 283 284 PetscFunctionBegin; 285 ierr = PetscNew(Mat_IS,&b);CHKERRQ(ierr); 286 A->data = (void*)b; 287 ierr = PetscMemzero(b,sizeof(Mat_IS));CHKERRQ(ierr); 288 ierr = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 289 A->factor = 0; 290 A->mapping = 0; 291 292 A->ops->mult = MatMult_IS; 293 A->ops->destroy = MatDestroy_IS; 294 A->ops->setlocaltoglobalmapping = MatSetLocalToGlobalMapping_IS; 295 A->ops->setvalueslocal = MatSetValuesLocal_IS; 296 A->ops->zerorowslocal = MatZeroRowsLocal_IS; 297 A->ops->assemblybegin = MatAssemblyBegin_IS; 298 A->ops->assemblyend = MatAssemblyEnd_IS; 299 A->ops->view = MatView_IS; 300 301 ierr = PetscSplitOwnership(A->comm,&A->m,&A->M);CHKERRQ(ierr); 302 ierr = PetscSplitOwnership(A->comm,&A->n,&A->N);CHKERRQ(ierr); 303 ierr = MPI_Scan(&A->m,&b->rend,1,MPI_INT,MPI_SUM,A->comm);CHKERRQ(ierr); 304 b->rstart = b->rend - A->m; 305 306 b->A = 0; 307 b->ctx = 0; 308 b->x = 0; 309 b->y = 0; 310 ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatISGetLocalMat_C","MatISGetLocalMat_IS",MatISGetLocalMat_IS);CHKERRQ(ierr); 311 312 PetscFunctionReturn(0); 313 } 314 EXTERN_C_END 315 316 317 318 319 320 321 322 323 324 325 326 327 328