17807a1faSBarry Smith 2273d9f13SBarry Smith #include "src/mat/matimpl.h" /*I "petscmat.h" I*/ 3325e03aeSBarry Smith #include "petscsys.h" 47807a1faSBarry Smith 54a2ae208SSatish Balay #undef __FUNCT__ 64a2ae208SSatish Balay #define __FUNCT__ "MatPublish_Base" 7*6849ba73SBarry Smith static PetscErrorCode MatPublish_Base(PetscObject obj) 835d8aa7fSBarry Smith { 935d8aa7fSBarry Smith #if defined(PETSC_HAVE_AMS) 1035d8aa7fSBarry Smith Mat mat = (Mat)obj; 11dfbe8321SBarry Smith PetscErrorCode ierr; 1235d8aa7fSBarry Smith #endif 1335d8aa7fSBarry Smith 1435d8aa7fSBarry Smith PetscFunctionBegin; 1535d8aa7fSBarry Smith #if defined(PETSC_HAVE_AMS) 1635d8aa7fSBarry Smith /* if it is already published then return */ 1735d8aa7fSBarry Smith if (mat->amem >=0) PetscFunctionReturn(0); 1835d8aa7fSBarry Smith 1935d8aa7fSBarry Smith ierr = PetscObjectPublishBaseBegin(obj);CHKERRQ(ierr); 2035d8aa7fSBarry Smith ierr = AMS_Memory_add_field((AMS_Memory)mat->amem,"globalsizes",&mat->M,2,AMS_INT,AMS_READ, 2135d8aa7fSBarry Smith AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr); 2235d8aa7fSBarry Smith ierr = AMS_Memory_add_field((AMS_Memory)mat->amem,"localsizes",&mat->m,2,AMS_INT,AMS_READ, 2335d8aa7fSBarry Smith AMS_DISTRIBUTED,AMS_REDUCT_UNDEF);CHKERRQ(ierr); 2435d8aa7fSBarry Smith ierr = PetscObjectPublishBaseEnd(obj);CHKERRQ(ierr); 2535d8aa7fSBarry Smith #endif 2635d8aa7fSBarry Smith 2735d8aa7fSBarry Smith PetscFunctionReturn(0); 2835d8aa7fSBarry Smith } 2935d8aa7fSBarry Smith 3035d8aa7fSBarry Smith 314a2ae208SSatish Balay #undef __FUNCT__ 324a2ae208SSatish Balay #define __FUNCT__ "MatCreate" 33325ab940SBarry Smith /*@C 3469dd0797SLois Curfman McInnes MatCreate - Creates a matrix where the type is determined 357e5f4302SBarry Smith from either a call to MatSetType() or from the options database 367e5f4302SBarry Smith with a call to MatSetFromOptions(). The default matrix type is 377e5f4302SBarry Smith AIJ, using the routines MatCreateSeqAIJ() or MatCreateMPIAIJ() 387e5f4302SBarry Smith if you do not set a type in the options database. If you never 397e5f4302SBarry Smith call MatSetType() or MatSetFromOptions() it will generate an 40f8ab6608SSatish Balay error when you try to use the matrix. 4183e1b59cSLois Curfman McInnes 42cb13003dSBarry Smith Collective on MPI_Comm 43cb13003dSBarry Smith 447807a1faSBarry Smith Input Parameters: 4582b900a8SBarry Smith + m - number of local rows (or PETSC_DECIDE) 4682b900a8SBarry Smith . n - number of local columns (or PETSC_DECIDE) 4782b900a8SBarry Smith . M - number of global rows (or PETSC_DETERMINE) 4882b900a8SBarry Smith . N - number of global columns (or PETSC_DETERMINE) 49cb13003dSBarry Smith - comm - MPI communicator 507807a1faSBarry Smith 517807a1faSBarry Smith Output Parameter: 52dc401e71SLois Curfman McInnes . A - the matrix 53e0b365e2SLois Curfman McInnes 54273d9f13SBarry Smith Options Database Keys: 55273d9f13SBarry Smith + -mat_type seqaij - AIJ type, uses MatCreateSeqAIJ() 56273d9f13SBarry Smith . -mat_type mpiaij - AIJ type, uses MatCreateMPIAIJ() 57273d9f13SBarry Smith . -mat_type seqbdiag - block diagonal type, uses MatCreateSeqBDiag() 58273d9f13SBarry Smith . -mat_type mpibdiag - block diagonal type, uses MatCreateMPIBDiag() 59273d9f13SBarry Smith . -mat_type mpirowbs - rowbs type, uses MatCreateMPIRowbs() 60273d9f13SBarry Smith . -mat_type seqdense - dense type, uses MatCreateSeqDense() 61273d9f13SBarry Smith . -mat_type mpidense - dense type, uses MatCreateMPIDense() 62273d9f13SBarry Smith . -mat_type seqbaij - block AIJ type, uses MatCreateSeqBAIJ() 63273d9f13SBarry Smith - -mat_type mpibaij - block AIJ type, uses MatCreateMPIBAIJ() 64e0b365e2SLois Curfman McInnes 6583e1b59cSLois Curfman McInnes Even More Options Database Keys: 6683e1b59cSLois Curfman McInnes See the manpages for particular formats (e.g., MatCreateSeqAIJ()) 6783e1b59cSLois Curfman McInnes for additional format-specific options. 68e0b365e2SLois Curfman McInnes 69bd9ce289SLois Curfman McInnes Notes: 70ec6e0d80SSatish Balay If PETSC_DECIDE is not used for the arguments 'm' and 'n', then the 71ec6e0d80SSatish Balay user must ensure that they are chosen to be compatible with the 72ec6e0d80SSatish Balay vectors. To do this, one first considers the matrix-vector product 73ec6e0d80SSatish Balay 'y = A x'. The 'm' that is used in the above routine must match the 74ec6e0d80SSatish Balay local size used in the vector creation routine VecCreateMPI() for 'y'. 75ec6e0d80SSatish Balay Likewise, the 'n' used must match that used as the local size in 76ec6e0d80SSatish Balay VecCreateMPI() for 'x'. 77ec6e0d80SSatish Balay 78273d9f13SBarry Smith Level: beginner 79273d9f13SBarry Smith 80a2fc510eSBarry Smith User manual sections: 81a2fc510eSBarry Smith + sec_matcreate 82a2fc510eSBarry Smith - chapter_matrices 83a2fc510eSBarry Smith 84273d9f13SBarry Smith .keywords: matrix, create 85273d9f13SBarry Smith 86273d9f13SBarry Smith .seealso: MatCreateSeqAIJ((), MatCreateMPIAIJ(), 87273d9f13SBarry Smith MatCreateSeqBDiag(),MatCreateMPIBDiag(), 88273d9f13SBarry Smith MatCreateSeqDense(), MatCreateMPIDense(), 89273d9f13SBarry Smith MatCreateMPIRowbs(), MatCreateSeqBAIJ(), MatCreateMPIBAIJ(), 90273d9f13SBarry Smith MatCreateSeqSBAIJ(), MatCreateMPISBAIJ(), 91273d9f13SBarry Smith MatConvert() 92273d9f13SBarry Smith @*/ 93dfbe8321SBarry Smith PetscErrorCode MatCreate(MPI_Comm comm,int m,int n,int M,int N,Mat *A) 94273d9f13SBarry Smith { 95273d9f13SBarry Smith Mat B; 968ba1e511SMatthew Knepley #ifndef PETSC_USE_DYNAMIC_LIBRARIES 97dfbe8321SBarry Smith PetscErrorCode ierr; 988ba1e511SMatthew Knepley #endif 99273d9f13SBarry Smith 100273d9f13SBarry Smith PetscFunctionBegin; 1014482741eSBarry Smith PetscValidPointer(A,6); 1028ba1e511SMatthew Knepley *A = PETSC_NULL; 1038ba1e511SMatthew Knepley #ifndef PETSC_USE_DYNAMIC_LIBRARIES 1048ba1e511SMatthew Knepley ierr = MatInitializePackage(PETSC_NULL);CHKERRQ(ierr); 1058ba1e511SMatthew Knepley #endif 1068ba1e511SMatthew Knepley 107273d9f13SBarry Smith PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_COOKIE,0,"Mat",comm,MatDestroy,MatView); 108b0a32e0cSBarry Smith PetscLogObjectCreate(B); 109273d9f13SBarry Smith 110273d9f13SBarry Smith B->m = m; 111273d9f13SBarry Smith B->n = n; 112273d9f13SBarry Smith B->M = M; 113273d9f13SBarry Smith B->N = N; 114273d9f13SBarry Smith 115273d9f13SBarry Smith B->preallocated = PETSC_FALSE; 11635d8aa7fSBarry Smith B->bops->publish = MatPublish_Base; 117273d9f13SBarry Smith *A = B; 118273d9f13SBarry Smith PetscFunctionReturn(0); 119273d9f13SBarry Smith } 120273d9f13SBarry Smith 1214a2ae208SSatish Balay #undef __FUNCT__ 1224a2ae208SSatish Balay #define __FUNCT__ "MatSetFromOptions" 123273d9f13SBarry Smith /*@C 124273d9f13SBarry Smith MatSetFromOptions - Creates a matrix where the type is determined 125273d9f13SBarry Smith from the options database. Generates a parallel MPI matrix if the 126273d9f13SBarry Smith communicator has more than one processor. The default matrix type is 1277e5f4302SBarry Smith AIJ, using the routines MatCreateSeqAIJ() and MatCreateMPIAIJ() if 1287e5f4302SBarry Smith you do not select a type in the options database. 129273d9f13SBarry Smith 130273d9f13SBarry Smith Collective on Mat 131273d9f13SBarry Smith 132273d9f13SBarry Smith Input Parameter: 133273d9f13SBarry Smith . A - the matrix 134273d9f13SBarry Smith 135273d9f13SBarry Smith Options Database Keys: 136273d9f13SBarry Smith + -mat_type seqaij - AIJ type, uses MatCreateSeqAIJ() 137273d9f13SBarry Smith . -mat_type mpiaij - AIJ type, uses MatCreateMPIAIJ() 138273d9f13SBarry Smith . -mat_type seqbdiag - block diagonal type, uses MatCreateSeqBDiag() 139273d9f13SBarry Smith . -mat_type mpibdiag - block diagonal type, uses MatCreateMPIBDiag() 140273d9f13SBarry Smith . -mat_type mpirowbs - rowbs type, uses MatCreateMPIRowbs() 141273d9f13SBarry Smith . -mat_type seqdense - dense type, uses MatCreateSeqDense() 142273d9f13SBarry Smith . -mat_type mpidense - dense type, uses MatCreateMPIDense() 143273d9f13SBarry Smith . -mat_type seqbaij - block AIJ type, uses MatCreateSeqBAIJ() 144273d9f13SBarry Smith - -mat_type mpibaij - block AIJ type, uses MatCreateMPIBAIJ() 145273d9f13SBarry Smith 146273d9f13SBarry Smith Even More Options Database Keys: 147273d9f13SBarry Smith See the manpages for particular formats (e.g., MatCreateSeqAIJ()) 148273d9f13SBarry Smith for additional format-specific options. 149bd9ce289SLois Curfman McInnes 1501d69843bSLois Curfman McInnes Level: beginner 1511d69843bSLois Curfman McInnes 152dc401e71SLois Curfman McInnes .keywords: matrix, create 153e0b365e2SLois Curfman McInnes 154fafbff53SBarry Smith .seealso: MatCreateSeqAIJ((), MatCreateMPIAIJ(), 155fafbff53SBarry Smith MatCreateSeqBDiag(),MatCreateMPIBDiag(), 15639ddd567SLois Curfman McInnes MatCreateSeqDense(), MatCreateMPIDense(), 157a209d233SLois Curfman McInnes MatCreateMPIRowbs(), MatCreateSeqBAIJ(), MatCreateMPIBAIJ(), 158a209d233SLois Curfman McInnes MatCreateSeqSBAIJ(), MatCreateMPISBAIJ(), 159273d9f13SBarry Smith MatConvert() 1607807a1faSBarry Smith @*/ 161dfbe8321SBarry Smith PetscErrorCode MatSetFromOptions(Mat B) 1627807a1faSBarry Smith { 163dfbe8321SBarry Smith PetscErrorCode ierr; 164dfbe8321SBarry Smith int size; 165273d9f13SBarry Smith char mtype[256]; 166273d9f13SBarry Smith PetscTruth flg; 167dbb450caSBarry Smith 1683a40ed3dSBarry Smith PetscFunctionBegin; 169b0a32e0cSBarry Smith ierr = PetscOptionsGetString(B->prefix,"-mat_type",mtype,256,&flg);CHKERRQ(ierr); 170273d9f13SBarry Smith if (flg) { 171273d9f13SBarry Smith ierr = MatSetType(B,mtype);CHKERRQ(ierr); 172273d9f13SBarry Smith } 173273d9f13SBarry Smith if (!B->type_name) { 174273d9f13SBarry Smith ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr); 1751d9a38edSKris Buschelman ierr = MatSetType(B,MATAIJ);CHKERRQ(ierr); 176dfa27b74SSatish Balay } 177c7cd70f7SSatish Balay #if defined(__cplusplus) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE) && defined(PETSC_HAVE_CXX_NAMESPACE) 17824a595ddSBarry Smith ierr = MatESISetFromOptions(B);CHKERRQ(ierr); 17924a595ddSBarry Smith #endif 1803a40ed3dSBarry Smith PetscFunctionReturn(0); 1817807a1faSBarry Smith } 1827807a1faSBarry Smith 1834a2ae208SSatish Balay #undef __FUNCT__ 1844a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation" 185273d9f13SBarry Smith /*@C 186273d9f13SBarry Smith MatSetUpPreallocation 187dae03382SLois Curfman McInnes 188273d9f13SBarry Smith Collective on Mat 189dae03382SLois Curfman McInnes 190273d9f13SBarry Smith Input Parameter: 191273d9f13SBarry Smith . A - the matrix 192d5d45c9bSBarry Smith 193273d9f13SBarry Smith Level: beginner 194d5d45c9bSBarry Smith 195273d9f13SBarry Smith .keywords: matrix, create 196273d9f13SBarry Smith 197273d9f13SBarry Smith .seealso: MatCreateSeqAIJ((), MatCreateMPIAIJ(), 198273d9f13SBarry Smith MatCreateSeqBDiag(),MatCreateMPIBDiag(), 199273d9f13SBarry Smith MatCreateSeqDense(), MatCreateMPIDense(), 200273d9f13SBarry Smith MatCreateMPIRowbs(), MatCreateSeqBAIJ(), MatCreateMPIBAIJ(), 201273d9f13SBarry Smith MatCreateSeqSBAIJ(), MatCreateMPISBAIJ(), 202273d9f13SBarry Smith MatConvert() 203273d9f13SBarry Smith @*/ 204dfbe8321SBarry Smith PetscErrorCode MatSetUpPreallocation(Mat B) 205273d9f13SBarry Smith { 206dfbe8321SBarry Smith PetscErrorCode ierr; 207273d9f13SBarry Smith 208273d9f13SBarry Smith PetscFunctionBegin; 209273d9f13SBarry Smith if (B->ops->setuppreallocation) { 2108ba1e511SMatthew Knepley PetscLogInfo(B,"MatSetUpPreallocation: Warning not preallocating matrix storage"); 211273d9f13SBarry Smith ierr = (*B->ops->setuppreallocation)(B);CHKERRQ(ierr); 212273d9f13SBarry Smith B->ops->setuppreallocation = 0; 213273d9f13SBarry Smith B->preallocated = PETSC_TRUE; 214273d9f13SBarry Smith } 215273d9f13SBarry Smith PetscFunctionReturn(0); 216273d9f13SBarry Smith } 217273d9f13SBarry Smith 218273d9f13SBarry Smith /* 219273d9f13SBarry Smith Copies from Cs header to A 220273d9f13SBarry Smith */ 2214a2ae208SSatish Balay #undef __FUNCT__ 2224a2ae208SSatish Balay #define __FUNCT__ "MatHeaderCopy" 223dfbe8321SBarry Smith PetscErrorCode MatHeaderCopy(Mat A,Mat C) 224273d9f13SBarry Smith { 225dfbe8321SBarry Smith PetscErrorCode ierr; 226dfbe8321SBarry Smith int refct; 227273d9f13SBarry Smith PetscOps *Abops; 228273d9f13SBarry Smith MatOps Aops; 229273d9f13SBarry Smith char *mtype,*mname; 230273d9f13SBarry Smith 231273d9f13SBarry Smith PetscFunctionBegin; 232273d9f13SBarry Smith /* free all the interior data structures from mat */ 233273d9f13SBarry Smith ierr = (*A->ops->destroy)(A);CHKERRQ(ierr); 234273d9f13SBarry Smith 2358a124369SBarry Smith ierr = PetscMapDestroy(A->rmap);CHKERRQ(ierr); 2368a124369SBarry Smith ierr = PetscMapDestroy(A->cmap);CHKERRQ(ierr); 237273d9f13SBarry Smith 238273d9f13SBarry Smith /* save the parts of A we need */ 239273d9f13SBarry Smith Abops = A->bops; 240273d9f13SBarry Smith Aops = A->ops; 241273d9f13SBarry Smith refct = A->refct; 242273d9f13SBarry Smith mtype = A->type_name; 243273d9f13SBarry Smith mname = A->name; 244273d9f13SBarry Smith 245273d9f13SBarry Smith /* copy C over to A */ 246273d9f13SBarry Smith ierr = PetscMemcpy(A,C,sizeof(struct _p_Mat));CHKERRQ(ierr); 247273d9f13SBarry Smith 248273d9f13SBarry Smith /* return the parts of A we saved */ 249273d9f13SBarry Smith A->bops = Abops; 250273d9f13SBarry Smith A->ops = Aops; 251273d9f13SBarry Smith A->qlist = 0; 252273d9f13SBarry Smith A->refct = refct; 253273d9f13SBarry Smith A->type_name = mtype; 254273d9f13SBarry Smith A->name = mname; 255273d9f13SBarry Smith 2560cc1bc1eSSatish Balay PetscLogObjectDestroy(C); 257273d9f13SBarry Smith PetscHeaderDestroy(C); 258273d9f13SBarry Smith PetscFunctionReturn(0); 259273d9f13SBarry Smith } 260