1 #define PETSCMAT_DLL 2 3 #include "private/matimpl.h" /*I "petscmat.h" I*/ 4 5 #if 0 6 #undef __FUNCT__ 7 #define __FUNCT__ "MatPublish_Base" 8 static PetscErrorCode MatPublish_Base(PetscObject obj) 9 { 10 PetscFunctionBegin; 11 PetscFunctionReturn(0); 12 } 13 #endif 14 15 #undef __FUNCT__ 16 #define __FUNCT__ "MatCreate" 17 /*@ 18 MatCreate - Creates a matrix where the type is determined 19 from either a call to MatSetType() or from the options database 20 with a call to MatSetFromOptions(). The default matrix type is 21 AIJ, using the routines MatCreateSeqAIJ() or MatCreateMPIAIJ() 22 if you do not set a type in the options database. If you never 23 call MatSetType() or MatSetFromOptions() it will generate an 24 error when you try to use the matrix. 25 26 Collective on MPI_Comm 27 28 Input Parameter: 29 . comm - MPI communicator 30 31 Output Parameter: 32 . A - the matrix 33 34 Options Database Keys: 35 + -mat_type seqaij - AIJ type, uses MatCreateSeqAIJ() 36 . -mat_type mpiaij - AIJ type, uses MatCreateMPIAIJ() 37 . -mat_type seqdense - dense type, uses MatCreateSeqDense() 38 . -mat_type mpidense - dense type, uses MatCreateMPIDense() 39 . -mat_type seqbaij - block AIJ type, uses MatCreateSeqBAIJ() 40 - -mat_type mpibaij - block AIJ type, uses MatCreateMPIBAIJ() 41 42 Even More Options Database Keys: 43 See the manpages for particular formats (e.g., MatCreateSeqAIJ()) 44 for additional format-specific options. 45 46 Notes: 47 48 Level: beginner 49 50 User manual sections: 51 + sec_matcreate 52 - chapter_matrices 53 54 .keywords: matrix, create 55 56 .seealso: MatCreateSeqAIJ(), MatCreateMPIAIJ(), 57 MatCreateSeqDense(), MatCreateMPIDense(), 58 MatCreateSeqBAIJ(), MatCreateMPIBAIJ(), 59 MatCreateSeqSBAIJ(), MatCreateMPISBAIJ(), 60 MatConvert() 61 @*/ 62 PetscErrorCode PETSCMAT_DLLEXPORT MatCreate(MPI_Comm comm,Mat *A) 63 { 64 Mat B; 65 PetscErrorCode ierr; 66 67 PetscFunctionBegin; 68 PetscValidPointer(A,2); 69 70 *A = PETSC_NULL; 71 #ifndef PETSC_USE_DYNAMIC_LIBRARIES 72 ierr = MatInitializePackage(PETSC_NULL);CHKERRQ(ierr); 73 #endif 74 75 ierr = PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_CLASSID,0,"Mat",comm,MatDestroy,MatView);CHKERRQ(ierr); 76 ierr = PetscLayoutCreate(comm,&B->rmap);CHKERRQ(ierr); 77 ierr = PetscLayoutCreate(comm,&B->cmap);CHKERRQ(ierr); 78 B->preallocated = PETSC_FALSE; 79 *A = B; 80 PetscFunctionReturn(0); 81 } 82 83 #undef __FUNCT__ 84 #define __FUNCT__ "MatSetSizes" 85 /*@ 86 MatSetSizes - Sets the local and global sizes, and checks to determine compatibility 87 88 Collective on Mat 89 90 Input Parameters: 91 + A - the matrix 92 . m - number of local rows (or PETSC_DECIDE) 93 . n - number of local columns (or PETSC_DECIDE) 94 . M - number of global rows (or PETSC_DETERMINE) 95 - N - number of global columns (or PETSC_DETERMINE) 96 97 Notes: 98 m (n) and M (N) cannot be both PETSC_DECIDE 99 If one processor calls this with M (N) of PETSC_DECIDE then all processors must, otherwise the program will hang. 100 101 If PETSC_DECIDE is not used for the arguments 'm' and 'n', then the 102 user must ensure that they are chosen to be compatible with the 103 vectors. To do this, one first considers the matrix-vector product 104 'y = A x'. The 'm' that is used in the above routine must match the 105 local size used in the vector creation routine VecCreateMPI() for 'y'. 106 Likewise, the 'n' used must match that used as the local size in 107 VecCreateMPI() for 'x'. 108 109 Level: beginner 110 111 .seealso: MatGetSize(), PetscSplitOwnership() 112 @*/ 113 PetscErrorCode PETSCMAT_DLLEXPORT MatSetSizes(Mat A, PetscInt m, PetscInt n, PetscInt M, PetscInt N) 114 { 115 PetscErrorCode ierr; 116 117 PetscFunctionBegin; 118 PetscValidHeaderSpecific(A,MAT_CLASSID,1); 119 if (M > 0 && m > M) SETERRQ2(PETSC_ERR_ARG_INCOMP,"Local column size %D cannot be larger than global column size %D",m,M); 120 if (N > 0 && n > N) SETERRQ2(PETSC_ERR_ARG_INCOMP,"Local row size %D cannot be larger than global row size %D",n,N); 121 if (A->ops->setsizes) { 122 /* Since this will not be set until the type has been set, this will NOT be called on the initial 123 call of MatSetSizes() (which must be called BEFORE MatSetType() */ 124 ierr = (*A->ops->setsizes)(A,m,n,M,N);CHKERRQ(ierr); 125 } else { 126 if ((A->rmap->n >= 0 || A->rmap->N >= 0) && (A->rmap->n != m || A->rmap->N != M)) SETERRQ4(PETSC_ERR_SUP,"Cannot change/reset row sizes to %D local %D global after previously setting them to %D local %D global",m,M,A->rmap->n,A->rmap->N); 127 if ((A->cmap->n >= 0 || A->cmap->N >= 0) && (A->cmap->n != n || A->cmap->N != N)) SETERRQ4(PETSC_ERR_SUP,"Cannot change/reset column sizes to %D local %D global after previously setting them to %D local %D global",n,N,A->cmap->n,A->cmap->N); 128 } 129 A->rmap->n = m; 130 A->cmap->n = n; 131 A->rmap->N = M; 132 A->cmap->N = N; 133 if (A->ops->create) { 134 ierr = (*A->ops->create)(A);CHKERRQ(ierr); 135 A->ops->create = 0; 136 } 137 138 PetscFunctionReturn(0); 139 } 140 141 #undef __FUNCT__ 142 #define __FUNCT__ "MatSetFromOptions" 143 /*@ 144 MatSetFromOptions - Creates a matrix where the type is determined 145 from the options database. Generates a parallel MPI matrix if the 146 communicator has more than one processor. The default matrix type is 147 AIJ, using the routines MatCreateSeqAIJ() and MatCreateMPIAIJ() if 148 you do not select a type in the options database. 149 150 Collective on Mat 151 152 Input Parameter: 153 . A - the matrix 154 155 Options Database Keys: 156 + -mat_type seqaij - AIJ type, uses MatCreateSeqAIJ() 157 . -mat_type mpiaij - AIJ type, uses MatCreateMPIAIJ() 158 . -mat_type seqdense - dense type, uses MatCreateSeqDense() 159 . -mat_type mpidense - dense type, uses MatCreateMPIDense() 160 . -mat_type seqbaij - block AIJ type, uses MatCreateSeqBAIJ() 161 - -mat_type mpibaij - block AIJ type, uses MatCreateMPIBAIJ() 162 163 Even More Options Database Keys: 164 See the manpages for particular formats (e.g., MatCreateSeqAIJ()) 165 for additional format-specific options. 166 167 Level: beginner 168 169 .keywords: matrix, create 170 171 .seealso: MatCreateSeqAIJ((), MatCreateMPIAIJ(), 172 MatCreateSeqDense(), MatCreateMPIDense(), 173 MatCreateSeqBAIJ(), MatCreateMPIBAIJ(), 174 MatCreateSeqSBAIJ(), MatCreateMPISBAIJ(), 175 MatConvert() 176 @*/ 177 PetscErrorCode PETSCMAT_DLLEXPORT MatSetFromOptions(Mat B) 178 { 179 PetscErrorCode ierr; 180 const char *deft = MATAIJ; 181 char type[256]; 182 PetscTruth flg; 183 184 PetscFunctionBegin; 185 PetscValidHeaderSpecific(B,MAT_CLASSID,1); 186 187 ierr = PetscOptionsBegin(((PetscObject)B)->comm,((PetscObject)B)->prefix,"Matrix options","Mat");CHKERRQ(ierr); 188 ierr = PetscOptionsList("-mat_type","Matrix type","MatSetType",MatList,deft,type,256,&flg);CHKERRQ(ierr); 189 if (flg) { 190 ierr = MatSetType(B,type);CHKERRQ(ierr); 191 } else if (!((PetscObject)B)->type_name) { 192 ierr = MatSetType(B,deft);CHKERRQ(ierr); 193 } 194 195 if (B->ops->setfromoptions) { 196 ierr = (*B->ops->setfromoptions)(B);CHKERRQ(ierr); 197 } 198 199 ierr = PetscOptionsEnd();CHKERRQ(ierr); 200 201 PetscFunctionReturn(0); 202 } 203 204 #undef __FUNCT__ 205 #define __FUNCT__ "MatSetUpPreallocation" 206 /*@ 207 MatSetUpPreallocation 208 209 Collective on Mat 210 211 Input Parameter: 212 . A - the matrix 213 214 Level: beginner 215 216 .keywords: matrix, create 217 218 .seealso: MatCreateSeqAIJ((), MatCreateMPIAIJ(), 219 MatCreateSeqDense(), MatCreateMPIDense(), 220 MatCreateSeqBAIJ(), MatCreateMPIBAIJ(), 221 MatCreateSeqSBAIJ(), MatCreateMPISBAIJ(), 222 MatConvert() 223 @*/ 224 PetscErrorCode PETSCMAT_DLLEXPORT MatSetUpPreallocation(Mat B) 225 { 226 PetscErrorCode ierr; 227 228 PetscFunctionBegin; 229 if (!B->preallocated && B->ops->setuppreallocation) { 230 ierr = PetscInfo(B,"Warning not preallocating matrix storage\n");CHKERRQ(ierr); 231 ierr = (*B->ops->setuppreallocation)(B);CHKERRQ(ierr); 232 } 233 B->preallocated = PETSC_TRUE; 234 PetscFunctionReturn(0); 235 } 236 237 /* 238 Copies from Cs header to A 239 240 This is somewhat different from MatHeaderReplace() it would be nice to merge the code 241 */ 242 #undef __FUNCT__ 243 #define __FUNCT__ "MatHeaderCopy" 244 PetscErrorCode MatHeaderCopy(Mat A,Mat C) 245 { 246 PetscErrorCode ierr; 247 PetscInt refct; 248 PetscOps *Abops; 249 MatOps Aops; 250 char *mtype,*mname; 251 void *spptr; 252 253 PetscFunctionBegin; 254 /* save the parts of A we need */ 255 Abops = ((PetscObject)A)->bops; 256 Aops = A->ops; 257 refct = ((PetscObject)A)->refct; 258 mtype = ((PetscObject)A)->type_name; 259 mname = ((PetscObject)A)->name; 260 spptr = A->spptr; 261 262 /* zero these so the destroy below does not free them */ 263 ((PetscObject)A)->type_name = 0; 264 ((PetscObject)A)->name = 0; 265 266 /* free all the interior data structures from mat */ 267 ierr = (*A->ops->destroy)(A);CHKERRQ(ierr); 268 269 ierr = PetscFree(C->spptr);CHKERRQ(ierr); 270 271 ierr = PetscLayoutDestroy(A->rmap);CHKERRQ(ierr); 272 ierr = PetscLayoutDestroy(A->cmap);CHKERRQ(ierr); 273 ierr = PetscFListDestroy(&((PetscObject)A)->qlist);CHKERRQ(ierr); 274 ierr = PetscOListDestroy(((PetscObject)A)->olist);CHKERRQ(ierr); 275 276 /* copy C over to A */ 277 ierr = PetscMemcpy(A,C,sizeof(struct _p_Mat));CHKERRQ(ierr); 278 279 /* return the parts of A we saved */ 280 ((PetscObject)A)->bops = Abops; 281 A->ops = Aops; 282 ((PetscObject)A)->refct = refct; 283 ((PetscObject)A)->type_name = mtype; 284 ((PetscObject)A)->name = mname; 285 A->spptr = spptr; 286 287 /* since these two are copied into A we do not want them destroyed in C */ 288 ((PetscObject)C)->qlist = 0; 289 ((PetscObject)C)->olist = 0; 290 ierr = PetscHeaderDestroy(C);CHKERRQ(ierr); 291 PetscFunctionReturn(0); 292 } 293 /* 294 Replace A's header with that of C 295 This is essentially code moved from MatDestroy 296 297 This is somewhat different from MatHeaderCopy() it would be nice to merge the code 298 */ 299 #undef __FUNCT__ 300 #define __FUNCT__ "MatHeaderReplace" 301 PetscErrorCode MatHeaderReplace(Mat A,Mat C) 302 { 303 PetscErrorCode ierr; 304 305 PetscFunctionBegin; 306 if (A == C) PetscFunctionReturn(0); 307 308 /* free all the interior data structures from mat */ 309 ierr = (*A->ops->destroy)(A);CHKERRQ(ierr); 310 ierr = PetscHeaderDestroy_Private((PetscObject)A);CHKERRQ(ierr); 311 ierr = PetscFree(A->ops);CHKERRQ(ierr); 312 ierr = PetscLayoutDestroy(A->rmap);CHKERRQ(ierr); 313 ierr = PetscLayoutDestroy(A->cmap);CHKERRQ(ierr); 314 ierr = PetscFree(A->spptr);CHKERRQ(ierr); 315 316 /* copy C over to A */ 317 if (C) { 318 ierr = PetscMemcpy(A,C,sizeof(struct _p_Mat));CHKERRQ(ierr); 319 ierr = PetscLogObjectDestroy((PetscObject)C);CHKERRQ(ierr); 320 ierr = PetscFree(C);CHKERRQ(ierr); 321 } 322 PetscFunctionReturn(0); 323 } 324