1be1d678aSKris Buschelman #define PETSCMAT_DLL 27807a1faSBarry Smith 3273d9f13SBarry Smith #include "src/mat/matimpl.h" /*I "petscmat.h" I*/ 4325e03aeSBarry Smith #include "petscsys.h" 57807a1faSBarry Smith 64a2ae208SSatish Balay #undef __FUNCT__ 74a2ae208SSatish Balay #define __FUNCT__ "MatPublish_Base" 86849ba73SBarry Smith static PetscErrorCode MatPublish_Base(PetscObject obj) 935d8aa7fSBarry Smith { 1035d8aa7fSBarry Smith PetscFunctionBegin; 1135d8aa7fSBarry Smith PetscFunctionReturn(0); 1235d8aa7fSBarry Smith } 1335d8aa7fSBarry Smith 1435d8aa7fSBarry Smith 154a2ae208SSatish Balay #undef __FUNCT__ 164a2ae208SSatish Balay #define __FUNCT__ "MatCreate" 1705869f15SSatish Balay /*@ 1869dd0797SLois Curfman McInnes MatCreate - Creates a matrix where the type is determined 197e5f4302SBarry Smith from either a call to MatSetType() or from the options database 207e5f4302SBarry Smith with a call to MatSetFromOptions(). The default matrix type is 217e5f4302SBarry Smith AIJ, using the routines MatCreateSeqAIJ() or MatCreateMPIAIJ() 227e5f4302SBarry Smith if you do not set a type in the options database. If you never 237e5f4302SBarry Smith call MatSetType() or MatSetFromOptions() it will generate an 24f8ab6608SSatish Balay error when you try to use the matrix. 2583e1b59cSLois Curfman McInnes 26cb13003dSBarry Smith Collective on MPI_Comm 27cb13003dSBarry Smith 28f69a0ea3SMatthew Knepley Input Parameter: 29f69a0ea3SMatthew Knepley . comm - MPI communicator 307807a1faSBarry Smith 317807a1faSBarry Smith Output Parameter: 32dc401e71SLois Curfman McInnes . A - the matrix 33e0b365e2SLois Curfman McInnes 34273d9f13SBarry Smith Options Database Keys: 35273d9f13SBarry Smith + -mat_type seqaij - AIJ type, uses MatCreateSeqAIJ() 36273d9f13SBarry Smith . -mat_type mpiaij - AIJ type, uses MatCreateMPIAIJ() 37273d9f13SBarry Smith . -mat_type seqbdiag - block diagonal type, uses MatCreateSeqBDiag() 38273d9f13SBarry Smith . -mat_type mpibdiag - block diagonal type, uses MatCreateMPIBDiag() 39273d9f13SBarry Smith . -mat_type mpirowbs - rowbs type, uses MatCreateMPIRowbs() 40273d9f13SBarry Smith . -mat_type seqdense - dense type, uses MatCreateSeqDense() 41273d9f13SBarry Smith . -mat_type mpidense - dense type, uses MatCreateMPIDense() 42273d9f13SBarry Smith . -mat_type seqbaij - block AIJ type, uses MatCreateSeqBAIJ() 43273d9f13SBarry Smith - -mat_type mpibaij - block AIJ type, uses MatCreateMPIBAIJ() 44e0b365e2SLois Curfman McInnes 4583e1b59cSLois Curfman McInnes Even More Options Database Keys: 4683e1b59cSLois Curfman McInnes See the manpages for particular formats (e.g., MatCreateSeqAIJ()) 4783e1b59cSLois Curfman McInnes for additional format-specific options. 48e0b365e2SLois Curfman McInnes 49bd9ce289SLois Curfman McInnes Notes: 50ec6e0d80SSatish Balay 51273d9f13SBarry Smith Level: beginner 52273d9f13SBarry Smith 53a2fc510eSBarry Smith User manual sections: 54a2fc510eSBarry Smith + sec_matcreate 55a2fc510eSBarry Smith - chapter_matrices 56a2fc510eSBarry Smith 57273d9f13SBarry Smith .keywords: matrix, create 58273d9f13SBarry Smith 59273d9f13SBarry Smith .seealso: MatCreateSeqAIJ((), MatCreateMPIAIJ(), 60273d9f13SBarry Smith MatCreateSeqBDiag(),MatCreateMPIBDiag(), 61273d9f13SBarry Smith MatCreateSeqDense(), MatCreateMPIDense(), 62273d9f13SBarry Smith MatCreateMPIRowbs(), MatCreateSeqBAIJ(), MatCreateMPIBAIJ(), 63273d9f13SBarry Smith MatCreateSeqSBAIJ(), MatCreateMPISBAIJ(), 64273d9f13SBarry Smith MatConvert() 65273d9f13SBarry Smith @*/ 66f69a0ea3SMatthew Knepley PetscErrorCode PETSCMAT_DLLEXPORT MatCreate(MPI_Comm comm,Mat *A) 67273d9f13SBarry Smith { 68273d9f13SBarry Smith Mat B; 69dfbe8321SBarry Smith PetscErrorCode ierr; 70273d9f13SBarry Smith 71273d9f13SBarry Smith PetscFunctionBegin; 72f69a0ea3SMatthew Knepley PetscValidPointer(A,2); 7397f1f81fSBarry Smith 748ba1e511SMatthew Knepley *A = PETSC_NULL; 758ba1e511SMatthew Knepley #ifndef PETSC_USE_DYNAMIC_LIBRARIES 768ba1e511SMatthew Knepley ierr = MatInitializePackage(PETSC_NULL);CHKERRQ(ierr); 778ba1e511SMatthew Knepley #endif 788ba1e511SMatthew Knepley 7952e6d16bSBarry Smith ierr = PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_COOKIE,0,"Mat",comm,MatDestroy,MatView);CHKERRQ(ierr); 80899cda47SBarry Smith B->rmap.n = -1; 81899cda47SBarry Smith B->rmap.N = -1; 82899cda47SBarry Smith B->cmap.n = -1; 83899cda47SBarry Smith B->cmap.N = -1; 84899cda47SBarry Smith B->rmap.bs = 1; 85899cda47SBarry Smith B->cmap.bs = 1; 86273d9f13SBarry Smith B->preallocated = PETSC_FALSE; 8735d8aa7fSBarry Smith B->bops->publish = MatPublish_Base; 88273d9f13SBarry Smith *A = B; 89273d9f13SBarry Smith PetscFunctionReturn(0); 90273d9f13SBarry Smith } 91273d9f13SBarry Smith 924a2ae208SSatish Balay #undef __FUNCT__ 93f69a0ea3SMatthew Knepley #define __FUNCT__ "MatSetSizes" 94f69a0ea3SMatthew Knepley /*@ 95f69a0ea3SMatthew Knepley MatSetSizes - Sets the local and global sizes, and checks to determine compatibility 96f69a0ea3SMatthew Knepley 97f69a0ea3SMatthew Knepley Collective on Mat 98f69a0ea3SMatthew Knepley 99f69a0ea3SMatthew Knepley Input Parameters: 100f69a0ea3SMatthew Knepley + A - the matrix 101f69a0ea3SMatthew Knepley . m - number of local rows (or PETSC_DECIDE) 102f69a0ea3SMatthew Knepley . n - number of local columns (or PETSC_DECIDE) 103f69a0ea3SMatthew Knepley . M - number of global rows (or PETSC_DETERMINE) 104f69a0ea3SMatthew Knepley - N - number of global columns (or PETSC_DETERMINE) 105f69a0ea3SMatthew Knepley 106f69a0ea3SMatthew Knepley Notes: 107f69a0ea3SMatthew Knepley m (n) and M (N) cannot be both PETSC_DECIDE 108f69a0ea3SMatthew Knepley If one processor calls this with M (N) of PETSC_DECIDE then all processors must, otherwise the program will hang. 109f69a0ea3SMatthew Knepley 110f69a0ea3SMatthew Knepley If PETSC_DECIDE is not used for the arguments 'm' and 'n', then the 111f69a0ea3SMatthew Knepley user must ensure that they are chosen to be compatible with the 112f69a0ea3SMatthew Knepley vectors. To do this, one first considers the matrix-vector product 113f69a0ea3SMatthew Knepley 'y = A x'. The 'm' that is used in the above routine must match the 114f69a0ea3SMatthew Knepley local size used in the vector creation routine VecCreateMPI() for 'y'. 115f69a0ea3SMatthew Knepley Likewise, the 'n' used must match that used as the local size in 116f69a0ea3SMatthew Knepley VecCreateMPI() for 'x'. 117f69a0ea3SMatthew Knepley 118f69a0ea3SMatthew Knepley Level: beginner 119f69a0ea3SMatthew Knepley 120f69a0ea3SMatthew Knepley .seealso: MatGetSize(), PetscSplitOwnership() 121f69a0ea3SMatthew Knepley @*/ 122f69a0ea3SMatthew Knepley PetscErrorCode PETSCMAT_DLLEXPORT MatSetSizes(Mat A, PetscInt m, PetscInt n, PetscInt M, PetscInt N) 123f69a0ea3SMatthew Knepley { 124284134d9SBarry Smith PetscErrorCode ierr; 125284134d9SBarry Smith 126f69a0ea3SMatthew Knepley PetscFunctionBegin; 127f69a0ea3SMatthew Knepley PetscValidHeaderSpecific(A,MAT_COOKIE,1); 128f69a0ea3SMatthew Knepley if (M > 0 && m > M) SETERRQ2(PETSC_ERR_ARG_INCOMP,"Local column size %D cannot be larger than global column size %D",m,M); 129f69a0ea3SMatthew Knepley if (N > 0 && n > N) SETERRQ2(PETSC_ERR_ARG_INCOMP,"Local row size %D cannot be larger than global row size %D",n,N); 130284134d9SBarry Smith if (A->ops->setsizes) { 131284134d9SBarry Smith /* Since this will not be set until the type has been set, this will NOT be called on the initial 132284134d9SBarry Smith call of MatSetSizes() (which must be called BEFORE MatSetType() */ 133284134d9SBarry Smith ierr = (*A->ops->setsizes)(A,m,n,M,N);CHKERRQ(ierr); 134284134d9SBarry Smith } else { 135899cda47SBarry Smith 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); 136899cda47SBarry Smith 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); 137284134d9SBarry Smith } 138899cda47SBarry Smith A->rmap.n = m; 139899cda47SBarry Smith A->cmap.n = n; 140899cda47SBarry Smith A->rmap.N = M; 141899cda47SBarry Smith A->cmap.N = N; 142*9b5a856fSSatish Balay /* also set rmap/cmap.bs here? */ 143*9b5a856fSSatish Balay ierr = PetscMapInitialize(A->comm,&A->rmap);CHKERRQ(ierr); 144*9b5a856fSSatish Balay ierr = PetscMapInitialize(A->comm,&A->cmap);CHKERRQ(ierr); 145f69a0ea3SMatthew Knepley PetscFunctionReturn(0); 146f69a0ea3SMatthew Knepley } 147f69a0ea3SMatthew Knepley 148f69a0ea3SMatthew Knepley #undef __FUNCT__ 1494a2ae208SSatish Balay #define __FUNCT__ "MatSetFromOptions" 15005869f15SSatish Balay /*@ 151273d9f13SBarry Smith MatSetFromOptions - Creates a matrix where the type is determined 152273d9f13SBarry Smith from the options database. Generates a parallel MPI matrix if the 153273d9f13SBarry Smith communicator has more than one processor. The default matrix type is 1547e5f4302SBarry Smith AIJ, using the routines MatCreateSeqAIJ() and MatCreateMPIAIJ() if 1557e5f4302SBarry Smith you do not select a type in the options database. 156273d9f13SBarry Smith 157273d9f13SBarry Smith Collective on Mat 158273d9f13SBarry Smith 159273d9f13SBarry Smith Input Parameter: 160273d9f13SBarry Smith . A - the matrix 161273d9f13SBarry Smith 162273d9f13SBarry Smith Options Database Keys: 163273d9f13SBarry Smith + -mat_type seqaij - AIJ type, uses MatCreateSeqAIJ() 164273d9f13SBarry Smith . -mat_type mpiaij - AIJ type, uses MatCreateMPIAIJ() 165273d9f13SBarry Smith . -mat_type seqbdiag - block diagonal type, uses MatCreateSeqBDiag() 166273d9f13SBarry Smith . -mat_type mpibdiag - block diagonal type, uses MatCreateMPIBDiag() 167273d9f13SBarry Smith . -mat_type mpirowbs - rowbs type, uses MatCreateMPIRowbs() 168273d9f13SBarry Smith . -mat_type seqdense - dense type, uses MatCreateSeqDense() 169273d9f13SBarry Smith . -mat_type mpidense - dense type, uses MatCreateMPIDense() 170273d9f13SBarry Smith . -mat_type seqbaij - block AIJ type, uses MatCreateSeqBAIJ() 171273d9f13SBarry Smith - -mat_type mpibaij - block AIJ type, uses MatCreateMPIBAIJ() 172273d9f13SBarry Smith 173273d9f13SBarry Smith Even More Options Database Keys: 174273d9f13SBarry Smith See the manpages for particular formats (e.g., MatCreateSeqAIJ()) 175273d9f13SBarry Smith for additional format-specific options. 176bd9ce289SLois Curfman McInnes 1771d69843bSLois Curfman McInnes Level: beginner 1781d69843bSLois Curfman McInnes 179dc401e71SLois Curfman McInnes .keywords: matrix, create 180e0b365e2SLois Curfman McInnes 181fafbff53SBarry Smith .seealso: MatCreateSeqAIJ((), MatCreateMPIAIJ(), 182fafbff53SBarry Smith MatCreateSeqBDiag(),MatCreateMPIBDiag(), 18339ddd567SLois Curfman McInnes MatCreateSeqDense(), MatCreateMPIDense(), 184a209d233SLois Curfman McInnes MatCreateMPIRowbs(), MatCreateSeqBAIJ(), MatCreateMPIBAIJ(), 185a209d233SLois Curfman McInnes MatCreateSeqSBAIJ(), MatCreateMPISBAIJ(), 186273d9f13SBarry Smith MatConvert() 1877807a1faSBarry Smith @*/ 188be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSetFromOptions(Mat B) 1897807a1faSBarry Smith { 190dfbe8321SBarry Smith PetscErrorCode ierr; 191273d9f13SBarry Smith char mtype[256]; 192273d9f13SBarry Smith PetscTruth flg; 193dbb450caSBarry Smith 1943a40ed3dSBarry Smith PetscFunctionBegin; 195b0a32e0cSBarry Smith ierr = PetscOptionsGetString(B->prefix,"-mat_type",mtype,256,&flg);CHKERRQ(ierr); 196273d9f13SBarry Smith if (flg) { 197273d9f13SBarry Smith ierr = MatSetType(B,mtype);CHKERRQ(ierr); 198273d9f13SBarry Smith } 199273d9f13SBarry Smith if (!B->type_name) { 20030735b05SKris Buschelman ierr = MatSetType(B,MATAIJ);CHKERRQ(ierr); 201dfa27b74SSatish Balay } 2023a40ed3dSBarry Smith PetscFunctionReturn(0); 2037807a1faSBarry Smith } 2047807a1faSBarry Smith 2054a2ae208SSatish Balay #undef __FUNCT__ 2064a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation" 207273d9f13SBarry Smith /*@C 208273d9f13SBarry Smith MatSetUpPreallocation 209dae03382SLois Curfman McInnes 210273d9f13SBarry Smith Collective on Mat 211dae03382SLois Curfman McInnes 212273d9f13SBarry Smith Input Parameter: 213273d9f13SBarry Smith . A - the matrix 214d5d45c9bSBarry Smith 215273d9f13SBarry Smith Level: beginner 216d5d45c9bSBarry Smith 217273d9f13SBarry Smith .keywords: matrix, create 218273d9f13SBarry Smith 219273d9f13SBarry Smith .seealso: MatCreateSeqAIJ((), MatCreateMPIAIJ(), 220273d9f13SBarry Smith MatCreateSeqBDiag(),MatCreateMPIBDiag(), 221273d9f13SBarry Smith MatCreateSeqDense(), MatCreateMPIDense(), 222273d9f13SBarry Smith MatCreateMPIRowbs(), MatCreateSeqBAIJ(), MatCreateMPIBAIJ(), 223273d9f13SBarry Smith MatCreateSeqSBAIJ(), MatCreateMPISBAIJ(), 224273d9f13SBarry Smith MatConvert() 225273d9f13SBarry Smith @*/ 226be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSetUpPreallocation(Mat B) 227273d9f13SBarry Smith { 228dfbe8321SBarry Smith PetscErrorCode ierr; 229273d9f13SBarry Smith 230273d9f13SBarry Smith PetscFunctionBegin; 231273d9f13SBarry Smith if (B->ops->setuppreallocation) { 232ae15b995SBarry Smith ierr = PetscInfo(B,"Warning not preallocating matrix storage\n");CHKERRQ(ierr); 233273d9f13SBarry Smith ierr = (*B->ops->setuppreallocation)(B);CHKERRQ(ierr); 234273d9f13SBarry Smith B->ops->setuppreallocation = 0; 235273d9f13SBarry Smith } 236210f0121SBarry Smith B->preallocated = PETSC_TRUE; 237273d9f13SBarry Smith PetscFunctionReturn(0); 238273d9f13SBarry Smith } 239273d9f13SBarry Smith 240273d9f13SBarry Smith /* 241273d9f13SBarry Smith Copies from Cs header to A 242273d9f13SBarry Smith */ 2434a2ae208SSatish Balay #undef __FUNCT__ 2444a2ae208SSatish Balay #define __FUNCT__ "MatHeaderCopy" 245dfbe8321SBarry Smith PetscErrorCode MatHeaderCopy(Mat A,Mat C) 246273d9f13SBarry Smith { 247dfbe8321SBarry Smith PetscErrorCode ierr; 248d44834fbSBarry Smith PetscInt refct; 249273d9f13SBarry Smith PetscOps *Abops; 250273d9f13SBarry Smith MatOps Aops; 251273d9f13SBarry Smith char *mtype,*mname; 25230735b05SKris Buschelman void *spptr; 253273d9f13SBarry Smith 254273d9f13SBarry Smith PetscFunctionBegin; 255273d9f13SBarry Smith /* free all the interior data structures from mat */ 256273d9f13SBarry Smith ierr = (*A->ops->destroy)(A);CHKERRQ(ierr); 257273d9f13SBarry Smith 258273d9f13SBarry Smith /* save the parts of A we need */ 259273d9f13SBarry Smith Abops = A->bops; 260273d9f13SBarry Smith Aops = A->ops; 261273d9f13SBarry Smith refct = A->refct; 262273d9f13SBarry Smith mtype = A->type_name; 263273d9f13SBarry Smith mname = A->name; 26430735b05SKris Buschelman spptr = A->spptr; 26530735b05SKris Buschelman 26630735b05SKris Buschelman ierr = PetscFree(C->spptr);CHKERRQ(ierr); 26730735b05SKris Buschelman C->spptr = PETSC_NULL; 26805b42c5fSBarry Smith 269357abbc8SBarry Smith ierr = PetscFree(A->rmap.range);CHKERRQ(ierr); 270357abbc8SBarry Smith ierr = PetscFree(A->cmap.range);CHKERRQ(ierr); 271273d9f13SBarry Smith 272273d9f13SBarry Smith /* copy C over to A */ 273273d9f13SBarry Smith ierr = PetscMemcpy(A,C,sizeof(struct _p_Mat));CHKERRQ(ierr); 274273d9f13SBarry Smith 275273d9f13SBarry Smith /* return the parts of A we saved */ 276273d9f13SBarry Smith A->bops = Abops; 277273d9f13SBarry Smith A->ops = Aops; 278273d9f13SBarry Smith A->qlist = 0; 279273d9f13SBarry Smith A->refct = refct; 280273d9f13SBarry Smith A->type_name = mtype; 281273d9f13SBarry Smith A->name = mname; 28230735b05SKris Buschelman A->spptr = spptr; 283273d9f13SBarry Smith 284d38fa0fbSBarry Smith ierr = PetscHeaderDestroy(C);CHKERRQ(ierr); 285273d9f13SBarry Smith PetscFunctionReturn(0); 286273d9f13SBarry Smith } 2878ab5b326SKris Buschelman /* 2888ab5b326SKris Buschelman Replace A's header with that of C 2898ab5b326SKris Buschelman This is essentially code moved from MatDestroy 2908ab5b326SKris Buschelman */ 2918ab5b326SKris Buschelman #undef __FUNCT__ 2928ab5b326SKris Buschelman #define __FUNCT__ "MatHeaderReplace" 2938ab5b326SKris Buschelman PetscErrorCode MatHeaderReplace(Mat A,Mat C) 2948ab5b326SKris Buschelman { 2958ab5b326SKris Buschelman PetscErrorCode ierr; 2968ab5b326SKris Buschelman 2978ab5b326SKris Buschelman PetscFunctionBegin; 2988ab5b326SKris Buschelman /* free all the interior data structures from mat */ 2998ab5b326SKris Buschelman ierr = (*A->ops->destroy)(A);CHKERRQ(ierr); 300d2c95936SBarry Smith ierr = PetscHeaderDestroy_Private((PetscObject)A);CHKERRQ(ierr); 301357abbc8SBarry Smith ierr = PetscFree(A->rmap.range);CHKERRQ(ierr); 302357abbc8SBarry Smith ierr = PetscFree(A->cmap.range);CHKERRQ(ierr); 3038ab5b326SKris Buschelman 3048ab5b326SKris Buschelman /* copy C over to A */ 3058ab5b326SKris Buschelman if (C) { 3068ab5b326SKris Buschelman ierr = PetscMemcpy(A,C,sizeof(struct _p_Mat));CHKERRQ(ierr); 3078ab5b326SKris Buschelman ierr = PetscLogObjectDestroy((PetscObject)C);CHKERRQ(ierr); 3088ab5b326SKris Buschelman ierr = PetscFree(C);CHKERRQ(ierr); 3098ab5b326SKris Buschelman } 3108ab5b326SKris Buschelman PetscFunctionReturn(0); 3118ab5b326SKris Buschelman } 312