1*7e5f4302SBarry Smith /*$Id: gcreate.c,v 1.128 2001/04/25 19:22:18 bsmith Exp bsmith $*/ 27807a1faSBarry Smith 3e090d566SSatish Balay #include "petscsys.h" 4273d9f13SBarry Smith #include "src/mat/matimpl.h" /*I "petscmat.h" I*/ 57807a1faSBarry Smith 64a2ae208SSatish Balay #undef __FUNCT__ 74a2ae208SSatish Balay #define __FUNCT__ "MatPublish_Base" 835d8aa7fSBarry Smith static int MatPublish_Base(PetscObject obj) 935d8aa7fSBarry Smith { 1035d8aa7fSBarry Smith #if defined(PETSC_HAVE_AMS) 1135d8aa7fSBarry Smith Mat mat = (Mat)obj; 1235d8aa7fSBarry Smith int ierr; 1335d8aa7fSBarry Smith #endif 1435d8aa7fSBarry Smith 1535d8aa7fSBarry Smith PetscFunctionBegin; 1635d8aa7fSBarry Smith #if defined(PETSC_HAVE_AMS) 1735d8aa7fSBarry Smith /* if it is already published then return */ 1835d8aa7fSBarry Smith if (mat->amem >=0) PetscFunctionReturn(0); 1935d8aa7fSBarry Smith 2035d8aa7fSBarry Smith ierr = PetscObjectPublishBaseBegin(obj);CHKERRQ(ierr); 2135d8aa7fSBarry Smith ierr = AMS_Memory_add_field((AMS_Memory)mat->amem,"globalsizes",&mat->M,2,AMS_INT,AMS_READ, 2235d8aa7fSBarry Smith AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr); 2335d8aa7fSBarry Smith ierr = AMS_Memory_add_field((AMS_Memory)mat->amem,"localsizes",&mat->m,2,AMS_INT,AMS_READ, 2435d8aa7fSBarry Smith AMS_DISTRIBUTED,AMS_REDUCT_UNDEF);CHKERRQ(ierr); 2535d8aa7fSBarry Smith ierr = PetscObjectPublishBaseEnd(obj);CHKERRQ(ierr); 2635d8aa7fSBarry Smith #endif 2735d8aa7fSBarry Smith 2835d8aa7fSBarry Smith PetscFunctionReturn(0); 2935d8aa7fSBarry Smith } 3035d8aa7fSBarry Smith 3135d8aa7fSBarry Smith 324a2ae208SSatish Balay #undef __FUNCT__ 334a2ae208SSatish Balay #define __FUNCT__ "MatCreate" 34325ab940SBarry Smith /*@C 3569dd0797SLois Curfman McInnes MatCreate - Creates a matrix where the type is determined 36*7e5f4302SBarry Smith from either a call to MatSetType() or from the options database 37*7e5f4302SBarry Smith with a call to MatSetFromOptions(). The default matrix type is 38*7e5f4302SBarry Smith AIJ, using the routines MatCreateSeqAIJ() or MatCreateMPIAIJ() 39*7e5f4302SBarry Smith if you do not set a type in the options database. If you never 40*7e5f4302SBarry Smith call MatSetType() or MatSetFromOptions() it will generate an 41*7e5f4302SBarry Smith eorror when you try to use the matrix. 4283e1b59cSLois Curfman McInnes 43cb13003dSBarry Smith Collective on MPI_Comm 44cb13003dSBarry Smith 457807a1faSBarry Smith Input Parameters: 4682b900a8SBarry Smith + m - number of local rows (or PETSC_DECIDE) 4782b900a8SBarry Smith . n - number of local columns (or PETSC_DECIDE) 4882b900a8SBarry Smith . M - number of global rows (or PETSC_DETERMINE) 4982b900a8SBarry Smith . N - number of global columns (or PETSC_DETERMINE) 50cb13003dSBarry Smith - comm - MPI communicator 517807a1faSBarry Smith 527807a1faSBarry Smith Output Parameter: 53dc401e71SLois Curfman McInnes . A - the matrix 54e0b365e2SLois Curfman McInnes 55273d9f13SBarry Smith Options Database Keys: 56273d9f13SBarry Smith + -mat_type seqaij - AIJ type, uses MatCreateSeqAIJ() 57273d9f13SBarry Smith . -mat_type mpiaij - AIJ type, uses MatCreateMPIAIJ() 58273d9f13SBarry Smith . -mat_type seqbdiag - block diagonal type, uses MatCreateSeqBDiag() 59273d9f13SBarry Smith . -mat_type mpibdiag - block diagonal type, uses MatCreateMPIBDiag() 60273d9f13SBarry Smith . -mat_type mpirowbs - rowbs type, uses MatCreateMPIRowbs() 61273d9f13SBarry Smith . -mat_type seqdense - dense type, uses MatCreateSeqDense() 62273d9f13SBarry Smith . -mat_type mpidense - dense type, uses MatCreateMPIDense() 63273d9f13SBarry Smith . -mat_type seqbaij - block AIJ type, uses MatCreateSeqBAIJ() 64273d9f13SBarry Smith - -mat_type mpibaij - block AIJ type, uses MatCreateMPIBAIJ() 65e0b365e2SLois Curfman McInnes 6683e1b59cSLois Curfman McInnes Even More Options Database Keys: 6783e1b59cSLois Curfman McInnes See the manpages for particular formats (e.g., MatCreateSeqAIJ()) 6883e1b59cSLois Curfman McInnes for additional format-specific options. 69e0b365e2SLois Curfman McInnes 70bd9ce289SLois Curfman McInnes Notes: 71ec6e0d80SSatish Balay If PETSC_DECIDE is not used for the arguments 'm' and 'n', then the 72ec6e0d80SSatish Balay user must ensure that they are chosen to be compatible with the 73ec6e0d80SSatish Balay vectors. To do this, one first considers the matrix-vector product 74ec6e0d80SSatish Balay 'y = A x'. The 'm' that is used in the above routine must match the 75ec6e0d80SSatish Balay local size used in the vector creation routine VecCreateMPI() for 'y'. 76ec6e0d80SSatish Balay Likewise, the 'n' used must match that used as the local size in 77ec6e0d80SSatish Balay VecCreateMPI() for 'x'. 78ec6e0d80SSatish Balay 79273d9f13SBarry Smith Level: beginner 80273d9f13SBarry Smith 81273d9f13SBarry Smith .keywords: matrix, create 82273d9f13SBarry Smith 83273d9f13SBarry Smith .seealso: MatCreateSeqAIJ((), MatCreateMPIAIJ(), 84273d9f13SBarry Smith MatCreateSeqBDiag(),MatCreateMPIBDiag(), 85273d9f13SBarry Smith MatCreateSeqDense(), MatCreateMPIDense(), 86273d9f13SBarry Smith MatCreateMPIRowbs(), MatCreateSeqBAIJ(), MatCreateMPIBAIJ(), 87273d9f13SBarry Smith MatCreateSeqSBAIJ(), MatCreateMPISBAIJ(), 88273d9f13SBarry Smith MatConvert() 89273d9f13SBarry Smith @*/ 90273d9f13SBarry Smith int MatCreate(MPI_Comm comm,int m,int n,int M,int N,Mat *A) 91273d9f13SBarry Smith { 92273d9f13SBarry Smith Mat B; 93273d9f13SBarry Smith 94273d9f13SBarry Smith PetscFunctionBegin; 95273d9f13SBarry Smith PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_COOKIE,0,"Mat",comm,MatDestroy,MatView); 96b0a32e0cSBarry Smith PetscLogObjectCreate(B); 97273d9f13SBarry Smith 98273d9f13SBarry Smith B->m = m; 99273d9f13SBarry Smith B->n = n; 100273d9f13SBarry Smith B->M = M; 101273d9f13SBarry Smith B->N = N; 102273d9f13SBarry Smith 103273d9f13SBarry Smith B->preallocated = PETSC_FALSE; 10435d8aa7fSBarry Smith B->bops->publish = MatPublish_Base; 105273d9f13SBarry Smith *A = B; 106273d9f13SBarry Smith PetscFunctionReturn(0); 107273d9f13SBarry Smith } 108273d9f13SBarry Smith 1094a2ae208SSatish Balay #undef __FUNCT__ 1104a2ae208SSatish Balay #define __FUNCT__ "MatSetFromOptions" 111273d9f13SBarry Smith /*@C 112273d9f13SBarry Smith MatSetFromOptions - Creates a matrix where the type is determined 113273d9f13SBarry Smith from the options database. Generates a parallel MPI matrix if the 114273d9f13SBarry Smith communicator has more than one processor. The default matrix type is 115*7e5f4302SBarry Smith AIJ, using the routines MatCreateSeqAIJ() and MatCreateMPIAIJ() if 116*7e5f4302SBarry Smith you do not select a type in the options database. 117273d9f13SBarry Smith 118273d9f13SBarry Smith Collective on Mat 119273d9f13SBarry Smith 120273d9f13SBarry Smith Input Parameter: 121273d9f13SBarry Smith . A - the matrix 122273d9f13SBarry Smith 123273d9f13SBarry Smith Options Database Keys: 124273d9f13SBarry Smith + -mat_type seqaij - AIJ type, uses MatCreateSeqAIJ() 125273d9f13SBarry Smith . -mat_type mpiaij - AIJ type, uses MatCreateMPIAIJ() 126273d9f13SBarry Smith . -mat_type seqbdiag - block diagonal type, uses MatCreateSeqBDiag() 127273d9f13SBarry Smith . -mat_type mpibdiag - block diagonal type, uses MatCreateMPIBDiag() 128273d9f13SBarry Smith . -mat_type mpirowbs - rowbs type, uses MatCreateMPIRowbs() 129273d9f13SBarry Smith . -mat_type seqdense - dense type, uses MatCreateSeqDense() 130273d9f13SBarry Smith . -mat_type mpidense - dense type, uses MatCreateMPIDense() 131273d9f13SBarry Smith . -mat_type seqbaij - block AIJ type, uses MatCreateSeqBAIJ() 132273d9f13SBarry Smith - -mat_type mpibaij - block AIJ type, uses MatCreateMPIBAIJ() 133273d9f13SBarry Smith 134273d9f13SBarry Smith Even More Options Database Keys: 135273d9f13SBarry Smith See the manpages for particular formats (e.g., MatCreateSeqAIJ()) 136273d9f13SBarry Smith for additional format-specific options. 137bd9ce289SLois Curfman McInnes 1381d69843bSLois Curfman McInnes Level: beginner 1391d69843bSLois Curfman McInnes 140dc401e71SLois Curfman McInnes .keywords: matrix, create 141e0b365e2SLois Curfman McInnes 142fafbff53SBarry Smith .seealso: MatCreateSeqAIJ((), MatCreateMPIAIJ(), 143fafbff53SBarry Smith MatCreateSeqBDiag(),MatCreateMPIBDiag(), 14439ddd567SLois Curfman McInnes MatCreateSeqDense(), MatCreateMPIDense(), 145a209d233SLois Curfman McInnes MatCreateMPIRowbs(), MatCreateSeqBAIJ(), MatCreateMPIBAIJ(), 146a209d233SLois Curfman McInnes MatCreateSeqSBAIJ(), MatCreateMPISBAIJ(), 147273d9f13SBarry Smith MatConvert() 1487807a1faSBarry Smith @*/ 149273d9f13SBarry Smith int MatSetFromOptions(Mat B) 1507807a1faSBarry Smith { 151273d9f13SBarry Smith int ierr,size; 152273d9f13SBarry Smith char mtype[256]; 153273d9f13SBarry Smith PetscTruth flg; 154dbb450caSBarry Smith 1553a40ed3dSBarry Smith PetscFunctionBegin; 156b0a32e0cSBarry Smith ierr = PetscOptionsGetString(B->prefix,"-mat_type",mtype,256,&flg);CHKERRQ(ierr); 157273d9f13SBarry Smith if (flg) { 158273d9f13SBarry Smith ierr = MatSetType(B,mtype);CHKERRQ(ierr); 159273d9f13SBarry Smith } 160273d9f13SBarry Smith if (!B->type_name) { 161273d9f13SBarry Smith ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr); 162273d9f13SBarry Smith if (size == 1) { 163273d9f13SBarry Smith ierr = MatSetType(B,MATSEQAIJ);CHKERRQ(ierr); 164273d9f13SBarry Smith } else { 165273d9f13SBarry Smith ierr = MatSetType(B,MATMPIAIJ);CHKERRQ(ierr); 166273d9f13SBarry Smith } 167dfa27b74SSatish Balay } 1683a40ed3dSBarry Smith PetscFunctionReturn(0); 1697807a1faSBarry Smith } 1707807a1faSBarry Smith 1714a2ae208SSatish Balay #undef __FUNCT__ 1724a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation" 173273d9f13SBarry Smith /*@C 174273d9f13SBarry Smith MatSetUpPreallocation 175dae03382SLois Curfman McInnes 176273d9f13SBarry Smith Collective on Mat 177dae03382SLois Curfman McInnes 178273d9f13SBarry Smith Input Parameter: 179273d9f13SBarry Smith . A - the matrix 180d5d45c9bSBarry Smith 181273d9f13SBarry Smith Level: beginner 182d5d45c9bSBarry Smith 183273d9f13SBarry Smith .keywords: matrix, create 184273d9f13SBarry Smith 185273d9f13SBarry Smith .seealso: MatCreateSeqAIJ((), MatCreateMPIAIJ(), 186273d9f13SBarry Smith MatCreateSeqBDiag(),MatCreateMPIBDiag(), 187273d9f13SBarry Smith MatCreateSeqDense(), MatCreateMPIDense(), 188273d9f13SBarry Smith MatCreateMPIRowbs(), MatCreateSeqBAIJ(), MatCreateMPIBAIJ(), 189273d9f13SBarry Smith MatCreateSeqSBAIJ(), MatCreateMPISBAIJ(), 190273d9f13SBarry Smith MatConvert() 191273d9f13SBarry Smith @*/ 192273d9f13SBarry Smith int MatSetUpPreallocation(Mat B) 193273d9f13SBarry Smith { 194273d9f13SBarry Smith int ierr; 195273d9f13SBarry Smith 196273d9f13SBarry Smith PetscFunctionBegin; 197273d9f13SBarry Smith if (B->ops->setuppreallocation) { 198b0a32e0cSBarry Smith PetscLogInfo(B,"MatSetTpPreallocation: Warning not preallocating matrix storage"); 199273d9f13SBarry Smith ierr = (*B->ops->setuppreallocation)(B);CHKERRQ(ierr); 200273d9f13SBarry Smith B->ops->setuppreallocation = 0; 201273d9f13SBarry Smith B->preallocated = PETSC_TRUE; 202273d9f13SBarry Smith } 203273d9f13SBarry Smith PetscFunctionReturn(0); 204273d9f13SBarry Smith } 205273d9f13SBarry Smith 206273d9f13SBarry Smith /* 207273d9f13SBarry Smith Copies from Cs header to A 208273d9f13SBarry Smith */ 2094a2ae208SSatish Balay #undef __FUNCT__ 2104a2ae208SSatish Balay #define __FUNCT__ "MatHeaderCopy" 211273d9f13SBarry Smith int MatHeaderCopy(Mat A,Mat C) 212273d9f13SBarry Smith { 213273d9f13SBarry Smith int ierr,refct; 214273d9f13SBarry Smith PetscOps *Abops; 215273d9f13SBarry Smith MatOps Aops; 216273d9f13SBarry Smith char *mtype,*mname; 217273d9f13SBarry Smith 218273d9f13SBarry Smith PetscFunctionBegin; 219273d9f13SBarry Smith /* free all the interior data structures from mat */ 220273d9f13SBarry Smith ierr = (*A->ops->destroy)(A);CHKERRQ(ierr); 221273d9f13SBarry Smith 222273d9f13SBarry Smith ierr = MapDestroy(A->rmap);CHKERRQ(ierr); 223273d9f13SBarry Smith ierr = MapDestroy(A->cmap);CHKERRQ(ierr); 224273d9f13SBarry Smith 225273d9f13SBarry Smith /* save the parts of A we need */ 226273d9f13SBarry Smith Abops = A->bops; 227273d9f13SBarry Smith Aops = A->ops; 228273d9f13SBarry Smith refct = A->refct; 229273d9f13SBarry Smith mtype = A->type_name; 230273d9f13SBarry Smith mname = A->name; 231273d9f13SBarry Smith 232273d9f13SBarry Smith /* copy C over to A */ 233273d9f13SBarry Smith ierr = PetscMemcpy(A,C,sizeof(struct _p_Mat));CHKERRQ(ierr); 234273d9f13SBarry Smith 235273d9f13SBarry Smith /* return the parts of A we saved */ 236273d9f13SBarry Smith A->bops = Abops; 237273d9f13SBarry Smith A->ops = Aops; 238273d9f13SBarry Smith A->qlist = 0; 239273d9f13SBarry Smith A->refct = refct; 240273d9f13SBarry Smith A->type_name = mtype; 241273d9f13SBarry Smith A->name = mname; 242273d9f13SBarry Smith 243273d9f13SBarry Smith PetscHeaderDestroy(C); 244273d9f13SBarry Smith PetscFunctionReturn(0); 245273d9f13SBarry Smith } 246