173f4d377SMatthew Knepley /*$Id: gcreate.c,v 1.131 2001/07/20 21:22:13 bsmith Exp $*/ 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" 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 367e5f4302SBarry Smith from either a call to MatSetType() or from the options database 377e5f4302SBarry Smith with a call to MatSetFromOptions(). The default matrix type is 387e5f4302SBarry Smith AIJ, using the routines MatCreateSeqAIJ() or MatCreateMPIAIJ() 397e5f4302SBarry Smith if you do not set a type in the options database. If you never 407e5f4302SBarry Smith call MatSetType() or MatSetFromOptions() it will generate an 41f8ab6608SSatish Balay error 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 81*a2fc510eSBarry Smith User manual sections: 82*a2fc510eSBarry Smith + sec_matcreate 83*a2fc510eSBarry Smith - chapter_matrices 84*a2fc510eSBarry Smith 85273d9f13SBarry Smith .keywords: matrix, create 86273d9f13SBarry Smith 87273d9f13SBarry Smith .seealso: MatCreateSeqAIJ((), MatCreateMPIAIJ(), 88273d9f13SBarry Smith MatCreateSeqBDiag(),MatCreateMPIBDiag(), 89273d9f13SBarry Smith MatCreateSeqDense(), MatCreateMPIDense(), 90273d9f13SBarry Smith MatCreateMPIRowbs(), MatCreateSeqBAIJ(), MatCreateMPIBAIJ(), 91273d9f13SBarry Smith MatCreateSeqSBAIJ(), MatCreateMPISBAIJ(), 92273d9f13SBarry Smith MatConvert() 93273d9f13SBarry Smith @*/ 94273d9f13SBarry Smith int MatCreate(MPI_Comm comm,int m,int n,int M,int N,Mat *A) 95273d9f13SBarry Smith { 96273d9f13SBarry Smith Mat B; 978ba1e511SMatthew Knepley #ifndef PETSC_USE_DYNAMIC_LIBRARIES 988ba1e511SMatthew Knepley int ierr; 998ba1e511SMatthew Knepley #endif 100273d9f13SBarry Smith 101273d9f13SBarry Smith PetscFunctionBegin; 1028ba1e511SMatthew Knepley PetscValidPointer(A); 1038ba1e511SMatthew Knepley *A = PETSC_NULL; 1048ba1e511SMatthew Knepley #ifndef PETSC_USE_DYNAMIC_LIBRARIES 1058ba1e511SMatthew Knepley ierr = MatInitializePackage(PETSC_NULL); CHKERRQ(ierr); 1068ba1e511SMatthew Knepley #endif 1078ba1e511SMatthew Knepley 108273d9f13SBarry Smith PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_COOKIE,0,"Mat",comm,MatDestroy,MatView); 109b0a32e0cSBarry Smith PetscLogObjectCreate(B); 110273d9f13SBarry Smith 111273d9f13SBarry Smith B->m = m; 112273d9f13SBarry Smith B->n = n; 113273d9f13SBarry Smith B->M = M; 114273d9f13SBarry Smith B->N = N; 115273d9f13SBarry Smith 116273d9f13SBarry Smith B->preallocated = PETSC_FALSE; 11735d8aa7fSBarry Smith B->bops->publish = MatPublish_Base; 118273d9f13SBarry Smith *A = B; 119273d9f13SBarry Smith PetscFunctionReturn(0); 120273d9f13SBarry Smith } 121273d9f13SBarry Smith 1224a2ae208SSatish Balay #undef __FUNCT__ 1234a2ae208SSatish Balay #define __FUNCT__ "MatSetFromOptions" 124273d9f13SBarry Smith /*@C 125273d9f13SBarry Smith MatSetFromOptions - Creates a matrix where the type is determined 126273d9f13SBarry Smith from the options database. Generates a parallel MPI matrix if the 127273d9f13SBarry Smith communicator has more than one processor. The default matrix type is 1287e5f4302SBarry Smith AIJ, using the routines MatCreateSeqAIJ() and MatCreateMPIAIJ() if 1297e5f4302SBarry Smith you do not select a type in the options database. 130273d9f13SBarry Smith 131273d9f13SBarry Smith Collective on Mat 132273d9f13SBarry Smith 133273d9f13SBarry Smith Input Parameter: 134273d9f13SBarry Smith . A - the matrix 135273d9f13SBarry Smith 136273d9f13SBarry Smith Options Database Keys: 137273d9f13SBarry Smith + -mat_type seqaij - AIJ type, uses MatCreateSeqAIJ() 138273d9f13SBarry Smith . -mat_type mpiaij - AIJ type, uses MatCreateMPIAIJ() 139273d9f13SBarry Smith . -mat_type seqbdiag - block diagonal type, uses MatCreateSeqBDiag() 140273d9f13SBarry Smith . -mat_type mpibdiag - block diagonal type, uses MatCreateMPIBDiag() 141273d9f13SBarry Smith . -mat_type mpirowbs - rowbs type, uses MatCreateMPIRowbs() 142273d9f13SBarry Smith . -mat_type seqdense - dense type, uses MatCreateSeqDense() 143273d9f13SBarry Smith . -mat_type mpidense - dense type, uses MatCreateMPIDense() 144273d9f13SBarry Smith . -mat_type seqbaij - block AIJ type, uses MatCreateSeqBAIJ() 145273d9f13SBarry Smith - -mat_type mpibaij - block AIJ type, uses MatCreateMPIBAIJ() 146273d9f13SBarry Smith 147273d9f13SBarry Smith Even More Options Database Keys: 148273d9f13SBarry Smith See the manpages for particular formats (e.g., MatCreateSeqAIJ()) 149273d9f13SBarry Smith for additional format-specific options. 150bd9ce289SLois Curfman McInnes 1511d69843bSLois Curfman McInnes Level: beginner 1521d69843bSLois Curfman McInnes 153dc401e71SLois Curfman McInnes .keywords: matrix, create 154e0b365e2SLois Curfman McInnes 155fafbff53SBarry Smith .seealso: MatCreateSeqAIJ((), MatCreateMPIAIJ(), 156fafbff53SBarry Smith MatCreateSeqBDiag(),MatCreateMPIBDiag(), 15739ddd567SLois Curfman McInnes MatCreateSeqDense(), MatCreateMPIDense(), 158a209d233SLois Curfman McInnes MatCreateMPIRowbs(), MatCreateSeqBAIJ(), MatCreateMPIBAIJ(), 159a209d233SLois Curfman McInnes MatCreateSeqSBAIJ(), MatCreateMPISBAIJ(), 160273d9f13SBarry Smith MatConvert() 1617807a1faSBarry Smith @*/ 162273d9f13SBarry Smith int MatSetFromOptions(Mat B) 1637807a1faSBarry Smith { 164273d9f13SBarry Smith int ierr,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); 175273d9f13SBarry Smith if (size == 1) { 176273d9f13SBarry Smith ierr = MatSetType(B,MATSEQAIJ);CHKERRQ(ierr); 177273d9f13SBarry Smith } else { 178273d9f13SBarry Smith ierr = MatSetType(B,MATMPIAIJ);CHKERRQ(ierr); 179273d9f13SBarry Smith } 180dfa27b74SSatish Balay } 181c1ef491bSBarry Smith #if defined(__cplusplus) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_MATSINGLE) && defined(PETSC_HAVE_CXX_NAMESPACE) 18224a595ddSBarry Smith ierr = MatESISetFromOptions(B);CHKERRQ(ierr); 18324a595ddSBarry Smith #endif 1843a40ed3dSBarry Smith PetscFunctionReturn(0); 1857807a1faSBarry Smith } 1867807a1faSBarry Smith 1874a2ae208SSatish Balay #undef __FUNCT__ 1884a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation" 189273d9f13SBarry Smith /*@C 190273d9f13SBarry Smith MatSetUpPreallocation 191dae03382SLois Curfman McInnes 192273d9f13SBarry Smith Collective on Mat 193dae03382SLois Curfman McInnes 194273d9f13SBarry Smith Input Parameter: 195273d9f13SBarry Smith . A - the matrix 196d5d45c9bSBarry Smith 197273d9f13SBarry Smith Level: beginner 198d5d45c9bSBarry Smith 199273d9f13SBarry Smith .keywords: matrix, create 200273d9f13SBarry Smith 201273d9f13SBarry Smith .seealso: MatCreateSeqAIJ((), MatCreateMPIAIJ(), 202273d9f13SBarry Smith MatCreateSeqBDiag(),MatCreateMPIBDiag(), 203273d9f13SBarry Smith MatCreateSeqDense(), MatCreateMPIDense(), 204273d9f13SBarry Smith MatCreateMPIRowbs(), MatCreateSeqBAIJ(), MatCreateMPIBAIJ(), 205273d9f13SBarry Smith MatCreateSeqSBAIJ(), MatCreateMPISBAIJ(), 206273d9f13SBarry Smith MatConvert() 207273d9f13SBarry Smith @*/ 208273d9f13SBarry Smith int MatSetUpPreallocation(Mat B) 209273d9f13SBarry Smith { 210273d9f13SBarry Smith int ierr; 211273d9f13SBarry Smith 212273d9f13SBarry Smith PetscFunctionBegin; 213273d9f13SBarry Smith if (B->ops->setuppreallocation) { 2148ba1e511SMatthew Knepley PetscLogInfo(B,"MatSetUpPreallocation: Warning not preallocating matrix storage"); 215273d9f13SBarry Smith ierr = (*B->ops->setuppreallocation)(B);CHKERRQ(ierr); 216273d9f13SBarry Smith B->ops->setuppreallocation = 0; 217273d9f13SBarry Smith B->preallocated = PETSC_TRUE; 218273d9f13SBarry Smith } 219273d9f13SBarry Smith PetscFunctionReturn(0); 220273d9f13SBarry Smith } 221273d9f13SBarry Smith 222273d9f13SBarry Smith /* 223273d9f13SBarry Smith Copies from Cs header to A 224273d9f13SBarry Smith */ 2254a2ae208SSatish Balay #undef __FUNCT__ 2264a2ae208SSatish Balay #define __FUNCT__ "MatHeaderCopy" 227273d9f13SBarry Smith int MatHeaderCopy(Mat A,Mat C) 228273d9f13SBarry Smith { 229273d9f13SBarry Smith int ierr,refct; 230273d9f13SBarry Smith PetscOps *Abops; 231273d9f13SBarry Smith MatOps Aops; 232273d9f13SBarry Smith char *mtype,*mname; 233273d9f13SBarry Smith 234273d9f13SBarry Smith PetscFunctionBegin; 235273d9f13SBarry Smith /* free all the interior data structures from mat */ 236273d9f13SBarry Smith ierr = (*A->ops->destroy)(A);CHKERRQ(ierr); 237273d9f13SBarry Smith 2388a124369SBarry Smith ierr = PetscMapDestroy(A->rmap);CHKERRQ(ierr); 2398a124369SBarry Smith ierr = PetscMapDestroy(A->cmap);CHKERRQ(ierr); 240273d9f13SBarry Smith 241273d9f13SBarry Smith /* save the parts of A we need */ 242273d9f13SBarry Smith Abops = A->bops; 243273d9f13SBarry Smith Aops = A->ops; 244273d9f13SBarry Smith refct = A->refct; 245273d9f13SBarry Smith mtype = A->type_name; 246273d9f13SBarry Smith mname = A->name; 247273d9f13SBarry Smith 248273d9f13SBarry Smith /* copy C over to A */ 249273d9f13SBarry Smith ierr = PetscMemcpy(A,C,sizeof(struct _p_Mat));CHKERRQ(ierr); 250273d9f13SBarry Smith 251273d9f13SBarry Smith /* return the parts of A we saved */ 252273d9f13SBarry Smith A->bops = Abops; 253273d9f13SBarry Smith A->ops = Aops; 254273d9f13SBarry Smith A->qlist = 0; 255273d9f13SBarry Smith A->refct = refct; 256273d9f13SBarry Smith A->type_name = mtype; 257273d9f13SBarry Smith A->name = mname; 258273d9f13SBarry Smith 2590cc1bc1eSSatish Balay PetscLogObjectDestroy(C); 260273d9f13SBarry Smith PetscHeaderDestroy(C); 261273d9f13SBarry Smith PetscFunctionReturn(0); 262273d9f13SBarry Smith } 263