1*bbf0e169SBarry Smith 2*bbf0e169SBarry Smith #ifndef lint 3*bbf0e169SBarry Smith static char vcid[] = "$Id: matrix.c,v 1.197 1996/10/02 17:22:28 bsmith Exp $"; 4*bbf0e169SBarry Smith #endif 5*bbf0e169SBarry Smith 6*bbf0e169SBarry Smith /* 7*bbf0e169SBarry Smith This is where the abstract matrix operations are defined 8*bbf0e169SBarry Smith */ 9*bbf0e169SBarry Smith 10*bbf0e169SBarry Smith #include "petsc.h" 11*bbf0e169SBarry Smith #include "src/mat/matimpl.h" /*I "mat.h" I*/ 12*bbf0e169SBarry Smith #include "src/vec/vecimpl.h" 13*bbf0e169SBarry Smith #include "pinclude/pviewer.h" 14*bbf0e169SBarry Smith #include "draw.h" 15*bbf0e169SBarry Smith 16*bbf0e169SBarry Smith 17*bbf0e169SBarry Smith /*@C 18*bbf0e169SBarry Smith MatGetRow - Gets a row of a matrix. You MUST call MatRestoreRow() 19*bbf0e169SBarry Smith for each row that you get to ensure that your application does 20*bbf0e169SBarry Smith not bleed memory. 21*bbf0e169SBarry Smith 22*bbf0e169SBarry Smith Input Parameters: 23*bbf0e169SBarry Smith . mat - the matrix 24*bbf0e169SBarry Smith . row - the row to get 25*bbf0e169SBarry Smith 26*bbf0e169SBarry Smith Output Parameters: 27*bbf0e169SBarry Smith . ncols - the number of nonzeros in the row 28*bbf0e169SBarry Smith . cols - if nonzero, the column numbers 29*bbf0e169SBarry Smith . vals - if nonzero, the values 30*bbf0e169SBarry Smith 31*bbf0e169SBarry Smith Notes: 32*bbf0e169SBarry Smith This routine is provided for people who need to have direct access 33*bbf0e169SBarry Smith to the structure of a matrix. We hope that we provide enough 34*bbf0e169SBarry Smith high-level matrix routines that few users will need it. 35*bbf0e169SBarry Smith 36*bbf0e169SBarry Smith For better efficiency, set cols and/or vals to PETSC_NULL if you do 37*bbf0e169SBarry Smith not wish to extract these quantities. 38*bbf0e169SBarry Smith 39*bbf0e169SBarry Smith The user can only examine the values extracted with MatGetRow(); 40*bbf0e169SBarry Smith the values cannot be altered. To change the matrix entries, one 41*bbf0e169SBarry Smith must use MatSetValues(). 42*bbf0e169SBarry Smith 43*bbf0e169SBarry Smith Caution: 44*bbf0e169SBarry Smith Do not try to change the contents of the output arrays (cols and vals). 45*bbf0e169SBarry Smith In some cases, this may corrupt the matrix. 46*bbf0e169SBarry Smith 47*bbf0e169SBarry Smith .keywords: matrix, row, get, extract 48*bbf0e169SBarry Smith 49*bbf0e169SBarry Smith .seealso: MatRestoreRow(), MatSetValues() 50*bbf0e169SBarry Smith @*/ 51*bbf0e169SBarry Smith int MatGetRow(Mat mat,int row,int *ncols,int **cols,Scalar **vals) 52*bbf0e169SBarry Smith { 53*bbf0e169SBarry Smith int ierr; 54*bbf0e169SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 55*bbf0e169SBarry Smith PetscValidIntPointer(ncols); 56*bbf0e169SBarry Smith if (!mat->assembled) SETERRQ(1,"MatGetRow:Not for unassembled matrix"); 57*bbf0e169SBarry Smith if (mat->factor) SETERRQ(1,"MatGetRow:Not for factored matrix"); 58*bbf0e169SBarry Smith PLogEventBegin(MAT_GetRow,mat,0,0,0); 59*bbf0e169SBarry Smith ierr = (*mat->ops.getrow)(mat,row,ncols,cols,vals); CHKERRQ(ierr); 60*bbf0e169SBarry Smith PLogEventEnd(MAT_GetRow,mat,0,0,0); 61*bbf0e169SBarry Smith return 0; 62*bbf0e169SBarry Smith } 63*bbf0e169SBarry Smith 64*bbf0e169SBarry Smith /*@C 65*bbf0e169SBarry Smith MatRestoreRow - Frees any temporary space allocated by MatGetRow(). 66*bbf0e169SBarry Smith 67*bbf0e169SBarry Smith Input Parameters: 68*bbf0e169SBarry Smith . mat - the matrix 69*bbf0e169SBarry Smith . row - the row to get 70*bbf0e169SBarry Smith . ncols, cols - the number of nonzeros and their columns 71*bbf0e169SBarry Smith . vals - if nonzero the column values 72*bbf0e169SBarry Smith 73*bbf0e169SBarry Smith .keywords: matrix, row, restore 74*bbf0e169SBarry Smith 75*bbf0e169SBarry Smith .seealso: MatGetRow() 76*bbf0e169SBarry Smith @*/ 77*bbf0e169SBarry Smith int MatRestoreRow(Mat mat,int row,int *ncols,int **cols,Scalar **vals) 78*bbf0e169SBarry Smith { 79*bbf0e169SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 80*bbf0e169SBarry Smith PetscValidIntPointer(ncols); 81*bbf0e169SBarry Smith if (!mat->assembled) SETERRQ(1,"MatRestoreRow:Not for unassembled matrix"); 82*bbf0e169SBarry Smith if (!mat->ops.restorerow) return 0; 83*bbf0e169SBarry Smith return (*mat->ops.restorerow)(mat,row,ncols,cols,vals); 84*bbf0e169SBarry Smith } 85*bbf0e169SBarry Smith /*@ 86*bbf0e169SBarry Smith MatView - Visualizes a matrix object. 87*bbf0e169SBarry Smith 88*bbf0e169SBarry Smith Input Parameters: 89*bbf0e169SBarry Smith . mat - the matrix 90*bbf0e169SBarry Smith . ptr - visualization context 91*bbf0e169SBarry Smith 92*bbf0e169SBarry Smith Notes: 93*bbf0e169SBarry Smith The available visualization contexts include 94*bbf0e169SBarry Smith $ VIEWER_STDOUT_SELF - standard output (default) 95*bbf0e169SBarry Smith $ VIEWER_STDOUT_WORLD - synchronized standard 96*bbf0e169SBarry Smith $ output where only the first processor opens 97*bbf0e169SBarry Smith $ the file. All other processors send their 98*bbf0e169SBarry Smith $ data to the first processor to print. 99*bbf0e169SBarry Smith 100*bbf0e169SBarry Smith The user can open alternative vistualization contexts with 101*bbf0e169SBarry Smith $ ViewerFileOpenASCII() - output to a specified file 102*bbf0e169SBarry Smith $ ViewerFileOpenBinary() - output in binary to a 103*bbf0e169SBarry Smith $ specified file; corresponding input uses MatLoad() 104*bbf0e169SBarry Smith $ ViewerDrawOpenX() - output nonzero matrix structure to 105*bbf0e169SBarry Smith $ an X window display 106*bbf0e169SBarry Smith $ ViewerMatlabOpen() - output matrix to Matlab viewer. 107*bbf0e169SBarry Smith $ Currently only the sequential dense and AIJ 108*bbf0e169SBarry Smith $ matrix types support the Matlab viewer. 109*bbf0e169SBarry Smith 110*bbf0e169SBarry Smith The user can call ViewerSetFormat() to specify the output 111*bbf0e169SBarry Smith format of ASCII printed objects (when using VIEWER_STDOUT_SELF, 112*bbf0e169SBarry Smith VIEWER_STDOUT_WORLD and ViewerFileOpenASCII). Available formats include 113*bbf0e169SBarry Smith $ VIEWER_FORMAT_ASCII_DEFAULT - default, prints matrix contents 114*bbf0e169SBarry Smith $ VIEWER_FORMAT_ASCII_MATLAB - Matlab format 115*bbf0e169SBarry Smith $ VIEWER_FORMAT_ASCII_IMPL - implementation-specific format 116*bbf0e169SBarry Smith $ (which is in many cases the same as the default) 117*bbf0e169SBarry Smith $ VIEWER_FORMAT_ASCII_INFO - basic information about the matrix 118*bbf0e169SBarry Smith $ size and structure (not the matrix entries) 119*bbf0e169SBarry Smith $ VIEWER_FORMAT_ASCII_INFO_LONG - more detailed information about the 120*bbf0e169SBarry Smith $ matrix structure 121*bbf0e169SBarry Smith 122*bbf0e169SBarry Smith .keywords: matrix, view, visualize, output, print, write, draw 123*bbf0e169SBarry Smith 124*bbf0e169SBarry Smith .seealso: ViewerSetFormat(), ViewerFileOpenASCII(), ViewerDrawOpenX(), 125*bbf0e169SBarry Smith ViewerMatlabOpen(), ViewerFileOpenBinary(), MatLoad() 126*bbf0e169SBarry Smith @*/ 127*bbf0e169SBarry Smith int MatView(Mat mat,Viewer viewer) 128*bbf0e169SBarry Smith { 129*bbf0e169SBarry Smith int format, ierr, rows, cols; 130*bbf0e169SBarry Smith FILE *fd; 131*bbf0e169SBarry Smith char *cstr; 132*bbf0e169SBarry Smith ViewerType vtype; 133*bbf0e169SBarry Smith MPI_Comm comm = mat->comm; 134*bbf0e169SBarry Smith 135*bbf0e169SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 136*bbf0e169SBarry Smith if (!mat->assembled) SETERRQ(1,"MatView:Not for unassembled matrix"); 137*bbf0e169SBarry Smith 138*bbf0e169SBarry Smith if (!viewer) { 139*bbf0e169SBarry Smith viewer = VIEWER_STDOUT_SELF; 140*bbf0e169SBarry Smith } 141*bbf0e169SBarry Smith 142*bbf0e169SBarry Smith ierr = ViewerGetType(viewer,&vtype); 143*bbf0e169SBarry Smith if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) { 144*bbf0e169SBarry Smith ierr = ViewerGetFormat(viewer,&format); CHKERRQ(ierr); 145*bbf0e169SBarry Smith ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr); 146*bbf0e169SBarry Smith if (format == VIEWER_FORMAT_ASCII_INFO || format == VIEWER_FORMAT_ASCII_INFO_LONG) { 147*bbf0e169SBarry Smith PetscFPrintf(comm,fd,"Matrix Object:\n"); 148*bbf0e169SBarry Smith ierr = MatGetType(mat,PETSC_NULL,&cstr); CHKERRQ(ierr); 149*bbf0e169SBarry Smith ierr = MatGetSize(mat,&rows,&cols); CHKERRQ(ierr); 150*bbf0e169SBarry Smith PetscFPrintf(comm,fd," type=%s, rows=%d, cols=%d\n",cstr,rows,cols); 151*bbf0e169SBarry Smith if (mat->ops.getinfo) { 152*bbf0e169SBarry Smith MatInfo info; 153*bbf0e169SBarry Smith ierr = MatGetInfo(mat,MAT_GLOBAL_SUM,&info); CHKERRQ(ierr); 154*bbf0e169SBarry Smith PetscFPrintf(comm,fd," total: nonzeros=%d, allocated nonzeros=%d\n", 155*bbf0e169SBarry Smith (int)info.nz_used,(int)info.nz_allocated); 156*bbf0e169SBarry Smith } 157*bbf0e169SBarry Smith } 158*bbf0e169SBarry Smith } 159*bbf0e169SBarry Smith if (mat->view) {ierr = (*mat->view)((PetscObject)mat,viewer); CHKERRQ(ierr);} 160*bbf0e169SBarry Smith return 0; 161*bbf0e169SBarry Smith } 162*bbf0e169SBarry Smith 163*bbf0e169SBarry Smith /*@C 164*bbf0e169SBarry Smith MatDestroy - Frees space taken by a matrix. 165*bbf0e169SBarry Smith 166*bbf0e169SBarry Smith Input Parameter: 167*bbf0e169SBarry Smith . mat - the matrix 168*bbf0e169SBarry Smith 169*bbf0e169SBarry Smith .keywords: matrix, destroy 170*bbf0e169SBarry Smith @*/ 171*bbf0e169SBarry Smith int MatDestroy(Mat mat) 172*bbf0e169SBarry Smith { 173*bbf0e169SBarry Smith int ierr; 174*bbf0e169SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 175*bbf0e169SBarry Smith ierr = (*mat->destroy)((PetscObject)mat); CHKERRQ(ierr); 176*bbf0e169SBarry Smith return 0; 177*bbf0e169SBarry Smith } 178*bbf0e169SBarry Smith /*@ 179*bbf0e169SBarry Smith MatValid - Checks whether a matrix object is valid. 180*bbf0e169SBarry Smith 181*bbf0e169SBarry Smith Input Parameter: 182*bbf0e169SBarry Smith . m - the matrix to check 183*bbf0e169SBarry Smith 184*bbf0e169SBarry Smith Output Parameter: 185*bbf0e169SBarry Smith flg - flag indicating matrix status, either 186*bbf0e169SBarry Smith $ PETSC_TRUE if matrix is valid; 187*bbf0e169SBarry Smith $ PETSC_FALSE otherwise. 188*bbf0e169SBarry Smith 189*bbf0e169SBarry Smith .keywords: matrix, valid 190*bbf0e169SBarry Smith @*/ 191*bbf0e169SBarry Smith int MatValid(Mat m,PetscTruth *flg) 192*bbf0e169SBarry Smith { 193*bbf0e169SBarry Smith PetscValidIntPointer(flg); 194*bbf0e169SBarry Smith if (!m) *flg = PETSC_FALSE; 195*bbf0e169SBarry Smith else if (m->cookie != MAT_COOKIE) *flg = PETSC_FALSE; 196*bbf0e169SBarry Smith else *flg = PETSC_TRUE; 197*bbf0e169SBarry Smith return 0; 198*bbf0e169SBarry Smith } 199*bbf0e169SBarry Smith 200*bbf0e169SBarry Smith /*@ 201*bbf0e169SBarry Smith MatSetValues - Inserts or adds a block of values into a matrix. 202*bbf0e169SBarry Smith These values may be cached, so MatAssemblyBegin() and MatAssemblyEnd() 203*bbf0e169SBarry Smith MUST be called after all calls to MatSetValues() have been completed. 204*bbf0e169SBarry Smith 205*bbf0e169SBarry Smith Input Parameters: 206*bbf0e169SBarry Smith . mat - the matrix 207*bbf0e169SBarry Smith . v - a logically two-dimensional array of values 208*bbf0e169SBarry Smith . m, indexm - the number of rows and their global indices 209*bbf0e169SBarry Smith . n, indexn - the number of columns and their global indices 210*bbf0e169SBarry Smith . addv - either ADD_VALUES or INSERT_VALUES, where 211*bbf0e169SBarry Smith $ ADD_VALUES - adds values to any existing entries 212*bbf0e169SBarry Smith $ INSERT_VALUES - replaces existing entries with new values 213*bbf0e169SBarry Smith 214*bbf0e169SBarry Smith Notes: 215*bbf0e169SBarry Smith By default the values, v, are row-oriented and unsorted. 216*bbf0e169SBarry Smith See MatSetOptions() for other options. 217*bbf0e169SBarry Smith 218*bbf0e169SBarry Smith Calls to MatSetValues() with the INSERT_VALUES and ADD_VALUES 219*bbf0e169SBarry Smith options cannot be mixed without intervening calls to the assembly 220*bbf0e169SBarry Smith routines. 221*bbf0e169SBarry Smith 222*bbf0e169SBarry Smith .keywords: matrix, insert, add, set, values 223*bbf0e169SBarry Smith 224*bbf0e169SBarry Smith .seealso: MatSetOptions(), MatAssemblyBegin(), MatAssemblyEnd() 225*bbf0e169SBarry Smith @*/ 226*bbf0e169SBarry Smith int MatSetValues(Mat mat,int m,int *idxm,int n,int *idxn,Scalar *v,InsertMode addv) 227*bbf0e169SBarry Smith { 228*bbf0e169SBarry Smith int ierr; 229*bbf0e169SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 230*bbf0e169SBarry Smith if (!m || !n) return 0; /* no values to insert */ 231*bbf0e169SBarry Smith PetscValidIntPointer(idxm); 232*bbf0e169SBarry Smith PetscValidIntPointer(idxn); 233*bbf0e169SBarry Smith PetscValidScalarPointer(v); 234*bbf0e169SBarry Smith if (mat->factor) SETERRQ(1,"MatSetValues:Not for factored matrix"); 235*bbf0e169SBarry Smith 236*bbf0e169SBarry Smith if (mat->assembled) { 237*bbf0e169SBarry Smith mat->was_assembled = PETSC_TRUE; 238*bbf0e169SBarry Smith mat->assembled = PETSC_FALSE; 239*bbf0e169SBarry Smith } 240*bbf0e169SBarry Smith PLogEventBegin(MAT_SetValues,mat,0,0,0); 241*bbf0e169SBarry Smith ierr = (*mat->ops.setvalues)(mat,m,idxm,n,idxn,v,addv);CHKERRQ(ierr); 242*bbf0e169SBarry Smith PLogEventEnd(MAT_SetValues,mat,0,0,0); 243*bbf0e169SBarry Smith return 0; 244*bbf0e169SBarry Smith } 245*bbf0e169SBarry Smith 246*bbf0e169SBarry Smith /*@ 247*bbf0e169SBarry Smith MatGetValues - Gets a block of values from a matrix. 248*bbf0e169SBarry Smith 249*bbf0e169SBarry Smith Input Parameters: 250*bbf0e169SBarry Smith . mat - the matrix 251*bbf0e169SBarry Smith . v - a logically two-dimensional array for storing the values 252*bbf0e169SBarry Smith . m, indexm - the number of rows and their global indices 253*bbf0e169SBarry Smith . n, indexn - the number of columns and their global indices 254*bbf0e169SBarry Smith 255*bbf0e169SBarry Smith Notes: 256*bbf0e169SBarry Smith The user must allocate space (m*n Scalars) for the values, v. 257*bbf0e169SBarry Smith The values, v, are then returned in a row-oriented format, 258*bbf0e169SBarry Smith analogous to that used by default in MatSetValues(). 259*bbf0e169SBarry Smith 260*bbf0e169SBarry Smith .keywords: matrix, get, values 261*bbf0e169SBarry Smith 262*bbf0e169SBarry Smith .seealso: MatGetRow(), MatGetSubmatrices(), MatSetValues() 263*bbf0e169SBarry Smith @*/ 264*bbf0e169SBarry Smith int MatGetValues(Mat mat,int m,int *idxm,int n,int *idxn,Scalar *v) 265*bbf0e169SBarry Smith { 266*bbf0e169SBarry Smith int ierr; 267*bbf0e169SBarry Smith 268*bbf0e169SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 269*bbf0e169SBarry Smith PetscValidIntPointer(idxm); 270*bbf0e169SBarry Smith PetscValidIntPointer(idxn); 271*bbf0e169SBarry Smith PetscValidScalarPointer(v); 272*bbf0e169SBarry Smith if (!mat->assembled) SETERRQ(1,"MatGetValues:Not for unassembled matrix"); 273*bbf0e169SBarry Smith if (mat->factor) SETERRQ(1,"MatGetValues:Not for factored matrix"); 274*bbf0e169SBarry Smith 275*bbf0e169SBarry Smith PLogEventBegin(MAT_GetValues,mat,0,0,0); 276*bbf0e169SBarry Smith ierr = (*mat->ops.getvalues)(mat,m,idxm,n,idxn,v); CHKERRQ(ierr); 277*bbf0e169SBarry Smith PLogEventEnd(MAT_GetValues,mat,0,0,0); 278*bbf0e169SBarry Smith return 0; 279*bbf0e169SBarry Smith } 280*bbf0e169SBarry Smith 281*bbf0e169SBarry Smith /* --------------------------------------------------------*/ 282*bbf0e169SBarry Smith /*@ 283*bbf0e169SBarry Smith MatMult - Computes the matrix-vector product, y = Ax. 284*bbf0e169SBarry Smith 285*bbf0e169SBarry Smith Input Parameters: 286*bbf0e169SBarry Smith . mat - the matrix 287*bbf0e169SBarry Smith . x - the vector to be multilplied 288*bbf0e169SBarry Smith 289*bbf0e169SBarry Smith Output Parameters: 290*bbf0e169SBarry Smith . y - the result 291*bbf0e169SBarry Smith 292*bbf0e169SBarry Smith .keywords: matrix, multiply, matrix-vector product 293*bbf0e169SBarry Smith 294*bbf0e169SBarry Smith .seealso: MatMultTrans(), MatMultAdd(), MatMultTransAdd() 295*bbf0e169SBarry Smith @*/ 296*bbf0e169SBarry Smith int MatMult(Mat mat,Vec x,Vec y) 297*bbf0e169SBarry Smith { 298*bbf0e169SBarry Smith int ierr; 299*bbf0e169SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 300*bbf0e169SBarry Smith PetscValidHeaderSpecific(x,VEC_COOKIE);PetscValidHeaderSpecific(y,VEC_COOKIE); 301*bbf0e169SBarry Smith if (!mat->assembled) SETERRQ(1,"MatMult:Not for unassembled matrix"); 302*bbf0e169SBarry Smith if (mat->factor) SETERRQ(1,"MatMult:Not for factored matrix"); 303*bbf0e169SBarry Smith if (x == y) SETERRQ(1,"MatMult:x and y must be different vectors"); 304*bbf0e169SBarry Smith if (mat->N != x->N) SETERRQ(PETSC_ERR_SIZ,"MatMult:Mat mat,Vec x: global dim"); 305*bbf0e169SBarry Smith if (mat->M != y->N) SETERRQ(PETSC_ERR_SIZ,"MatMult:Mat mat,Vec y: global dim"); 306*bbf0e169SBarry Smith if (mat->m != y->n) SETERRQ(PETSC_ERR_SIZ,"MatMult:Mat mat,Vec y: local dim"); 307*bbf0e169SBarry Smith 308*bbf0e169SBarry Smith PLogEventBegin(MAT_Mult,mat,x,y,0); 309*bbf0e169SBarry Smith ierr = (*mat->ops.mult)(mat,x,y); CHKERRQ(ierr); 310*bbf0e169SBarry Smith PLogEventEnd(MAT_Mult,mat,x,y,0); 311*bbf0e169SBarry Smith 312*bbf0e169SBarry Smith return 0; 313*bbf0e169SBarry Smith } 314*bbf0e169SBarry Smith /*@ 315*bbf0e169SBarry Smith MatMultTrans - Computes matrix transpose times a vector. 316*bbf0e169SBarry Smith 317*bbf0e169SBarry Smith Input Parameters: 318*bbf0e169SBarry Smith . mat - the matrix 319*bbf0e169SBarry Smith . x - the vector to be multilplied 320*bbf0e169SBarry Smith 321*bbf0e169SBarry Smith Output Parameters: 322*bbf0e169SBarry Smith . y - the result 323*bbf0e169SBarry Smith 324*bbf0e169SBarry Smith .keywords: matrix, multiply, matrix-vector product, transpose 325*bbf0e169SBarry Smith 326*bbf0e169SBarry Smith .seealso: MatMult(), MatMultAdd(), MatMultTransAdd() 327*bbf0e169SBarry Smith @*/ 328*bbf0e169SBarry Smith int MatMultTrans(Mat mat,Vec x,Vec y) 329*bbf0e169SBarry Smith { 330*bbf0e169SBarry Smith int ierr; 331*bbf0e169SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 332*bbf0e169SBarry Smith PetscValidHeaderSpecific(x,VEC_COOKIE); PetscValidHeaderSpecific(y,VEC_COOKIE); 333*bbf0e169SBarry Smith if (!mat->assembled) SETERRQ(1,"MatMultTrans:Not for unassembled matrix"); 334*bbf0e169SBarry Smith if (mat->factor) SETERRQ(1,"MatMult:Not for factored matrix"); 335*bbf0e169SBarry Smith if (x == y) SETERRQ(1,"MatMultTrans:x and y must be different vectors"); 336*bbf0e169SBarry Smith if (mat->M != x->N) SETERRQ(PETSC_ERR_SIZ,"MatMultTrans:Mat mat,Vec x: global dim"); 337*bbf0e169SBarry Smith if (mat->N != y->N) SETERRQ(PETSC_ERR_SIZ,"MatMultTrans:Mat mat,Vec y: global dim"); 338*bbf0e169SBarry Smith PLogEventBegin(MAT_MultTrans,mat,x,y,0); 339*bbf0e169SBarry Smith ierr = (*mat->ops.multtrans)(mat,x,y); CHKERRQ(ierr); 340*bbf0e169SBarry Smith PLogEventEnd(MAT_MultTrans,mat,x,y,0); 341*bbf0e169SBarry Smith return 0; 342*bbf0e169SBarry Smith } 343*bbf0e169SBarry Smith /*@ 344*bbf0e169SBarry Smith MatMultAdd - Computes v3 = v2 + A * v1. 345*bbf0e169SBarry Smith 346*bbf0e169SBarry Smith Input Parameters: 347*bbf0e169SBarry Smith . mat - the matrix 348*bbf0e169SBarry Smith . v1, v2 - the vectors 349*bbf0e169SBarry Smith 350*bbf0e169SBarry Smith Output Parameters: 351*bbf0e169SBarry Smith . v3 - the result 352*bbf0e169SBarry Smith 353*bbf0e169SBarry Smith .keywords: matrix, multiply, matrix-vector product, add 354*bbf0e169SBarry Smith 355*bbf0e169SBarry Smith .seealso: MatMultTrans(), MatMult(), MatMultTransAdd() 356*bbf0e169SBarry Smith @*/ 357*bbf0e169SBarry Smith int MatMultAdd(Mat mat,Vec v1,Vec v2,Vec v3) 358*bbf0e169SBarry Smith { 359*bbf0e169SBarry Smith int ierr; 360*bbf0e169SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE);PetscValidHeaderSpecific(v1,VEC_COOKIE); 361*bbf0e169SBarry Smith PetscValidHeaderSpecific(v2,VEC_COOKIE); PetscValidHeaderSpecific(v3,VEC_COOKIE); 362*bbf0e169SBarry Smith if (!mat->assembled) SETERRQ(1,"MatMultAdd:Not for unassembled matrix"); 363*bbf0e169SBarry Smith if (mat->factor) SETERRQ(1,"MatMult:Not for factored matrix"); 364*bbf0e169SBarry Smith if (mat->N != v1->N) SETERRQ(PETSC_ERR_SIZ,"MatMultAdd:Mat mat,Vec v1: global dim"); 365*bbf0e169SBarry Smith if (mat->M != v2->N) SETERRQ(PETSC_ERR_SIZ,"MatMultAdd:Mat mat,Vec v2: global dim"); 366*bbf0e169SBarry Smith if (mat->M != v3->N) SETERRQ(PETSC_ERR_SIZ,"MatMultAdd:Mat mat,Vec v3: global dim"); 367*bbf0e169SBarry Smith if (mat->m != v3->n) SETERRQ(PETSC_ERR_SIZ,"MatMultAdd:Mat mat,Vec v3: local dim"); 368*bbf0e169SBarry Smith if (mat->m != v2->n) SETERRQ(PETSC_ERR_SIZ,"MatMultAdd:Mat mat,Vec v2: local dim"); 369*bbf0e169SBarry Smith 370*bbf0e169SBarry Smith PLogEventBegin(MAT_MultAdd,mat,v1,v2,v3); 371*bbf0e169SBarry Smith if (v1 == v3) SETERRQ(1,"MatMultAdd:v1 and v3 must be different vectors"); 372*bbf0e169SBarry Smith ierr = (*mat->ops.multadd)(mat,v1,v2,v3); CHKERRQ(ierr); 373*bbf0e169SBarry Smith PLogEventEnd(MAT_MultAdd,mat,v1,v2,v3); 374*bbf0e169SBarry Smith return 0; 375*bbf0e169SBarry Smith } 376*bbf0e169SBarry Smith /*@ 377*bbf0e169SBarry Smith MatMultTransAdd - Computes v3 = v2 + A' * v1. 378*bbf0e169SBarry Smith 379*bbf0e169SBarry Smith Input Parameters: 380*bbf0e169SBarry Smith . mat - the matrix 381*bbf0e169SBarry Smith . v1, v2 - the vectors 382*bbf0e169SBarry Smith 383*bbf0e169SBarry Smith Output Parameters: 384*bbf0e169SBarry Smith . v3 - the result 385*bbf0e169SBarry Smith 386*bbf0e169SBarry Smith .keywords: matrix, multiply, matrix-vector product, transpose, add 387*bbf0e169SBarry Smith 388*bbf0e169SBarry Smith .seealso: MatMultTrans(), MatMultAdd(), MatMult() 389*bbf0e169SBarry Smith @*/ 390*bbf0e169SBarry Smith int MatMultTransAdd(Mat mat,Vec v1,Vec v2,Vec v3) 391*bbf0e169SBarry Smith { 392*bbf0e169SBarry Smith int ierr; 393*bbf0e169SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE);PetscValidHeaderSpecific(v1,VEC_COOKIE); 394*bbf0e169SBarry Smith PetscValidHeaderSpecific(v2,VEC_COOKIE);PetscValidHeaderSpecific(v3,VEC_COOKIE); 395*bbf0e169SBarry Smith if (!mat->assembled) SETERRQ(1,"MatMultTransAdd:Not for unassembled matrix"); 396*bbf0e169SBarry Smith if (mat->factor) SETERRQ(1,"MatMult:Not for factored matrix"); 397*bbf0e169SBarry Smith if (!mat->ops.multtransadd) SETERRQ(PETSC_ERR_SUP,"MatMultTransAdd"); 398*bbf0e169SBarry Smith if (v1 == v3) SETERRQ(1,"MatMultTransAdd:v1 and v2 must be different vectors"); 399*bbf0e169SBarry Smith if (mat->M != v1->N) SETERRQ(PETSC_ERR_SIZ,"MatMultTransAdd:Mat mat,Vec v1: global dim"); 400*bbf0e169SBarry Smith if (mat->N != v2->N) SETERRQ(PETSC_ERR_SIZ,"MatMultTransAdd:Mat mat,Vec v2: global dim"); 401*bbf0e169SBarry Smith if (mat->N != v3->N) SETERRQ(PETSC_ERR_SIZ,"MatMultTransAdd:Mat mat,Vec v3: global dim"); 402*bbf0e169SBarry Smith 403*bbf0e169SBarry Smith PLogEventBegin(MAT_MultTransAdd,mat,v1,v2,v3); 404*bbf0e169SBarry Smith ierr = (*mat->ops.multtransadd)(mat,v1,v2,v3); CHKERRQ(ierr); 405*bbf0e169SBarry Smith PLogEventEnd(MAT_MultTransAdd,mat,v1,v2,v3); 406*bbf0e169SBarry Smith return 0; 407*bbf0e169SBarry Smith } 408*bbf0e169SBarry Smith /* ------------------------------------------------------------*/ 409*bbf0e169SBarry Smith /*@C 410*bbf0e169SBarry Smith MatGetInfo - Returns information about matrix storage (number of 411*bbf0e169SBarry Smith nonzeros, memory, etc.). 412*bbf0e169SBarry Smith 413*bbf0e169SBarry Smith Input Parameters: 414*bbf0e169SBarry Smith . mat - the matrix 415*bbf0e169SBarry Smith 416*bbf0e169SBarry Smith Output Parameters: 417*bbf0e169SBarry Smith . flag - flag indicating the type of parameters to be returned 418*bbf0e169SBarry Smith $ flag = MAT_LOCAL: local matrix 419*bbf0e169SBarry Smith $ flag = MAT_GLOBAL_MAX: maximum over all processors 420*bbf0e169SBarry Smith $ flag = MAT_GLOBAL_SUM: sum over all processors 421*bbf0e169SBarry Smith . info - matrix information context 422*bbf0e169SBarry Smith 423*bbf0e169SBarry Smith Notes: 424*bbf0e169SBarry Smith The MatInfo context contains a variety of matrix data, including 425*bbf0e169SBarry Smith number of nonzeros allocated and used, number of mallocs during 426*bbf0e169SBarry Smith matrix assembly, etc. Additional information for factored matrices 427*bbf0e169SBarry Smith is provided (such as the fill ratio, number of mallocs during 428*bbf0e169SBarry Smith factorization, etc.). Much of this info is printed to STDOUT 429*bbf0e169SBarry Smith when using the runtime options 430*bbf0e169SBarry Smith $ -log_info -mat_view_info 431*bbf0e169SBarry Smith 432*bbf0e169SBarry Smith Example for C/C++ Users: 433*bbf0e169SBarry Smith See the file $(PETSC_DIR)/include/mat.h for a complete list of 434*bbf0e169SBarry Smith data within the MatInfo context. For example, 435*bbf0e169SBarry Smith $ 436*bbf0e169SBarry Smith $ MatInfo *info; 437*bbf0e169SBarry Smith $ Mat A; 438*bbf0e169SBarry Smith $ double mal, nz_a, nz_u; 439*bbf0e169SBarry Smith $ 440*bbf0e169SBarry Smith $ MatGetInfo(A,MAT_LOCAL,&info); 441*bbf0e169SBarry Smith $ mal = info->mallocs; 442*bbf0e169SBarry Smith $ nz_a = info->nz_allocated; 443*bbf0e169SBarry Smith $ 444*bbf0e169SBarry Smith 445*bbf0e169SBarry Smith Example for Fortran Users: 446*bbf0e169SBarry Smith Fortran users should declare info as a double precision 447*bbf0e169SBarry Smith array of dimension MAT_INFO_SIZE, and then extract the parameters 448*bbf0e169SBarry Smith of interest. See the file $(PETSC_DIR)/include/FINCLUDE/mat.h 449*bbf0e169SBarry Smith a complete list of parameter names. 450*bbf0e169SBarry Smith $ 451*bbf0e169SBarry Smith $ double precision info(MAT_INFO_SIZE) 452*bbf0e169SBarry Smith $ double precision mal, nz_a 453*bbf0e169SBarry Smith $ Mat A 454*bbf0e169SBarry Smith $ integer ierr 455*bbf0e169SBarry Smith $ 456*bbf0e169SBarry Smith $ call MatGetInfo(A,MAT_LOCAL,info,ierr) 457*bbf0e169SBarry Smith $ mal = info(MAT_INFO_MALLOCS) 458*bbf0e169SBarry Smith $ nz_a = info(MAT_INFO_NZ_ALLOCATED) 459*bbf0e169SBarry Smith $ 460*bbf0e169SBarry Smith 461*bbf0e169SBarry Smith .keywords: matrix, get, info, storage, nonzeros, memory, fill 462*bbf0e169SBarry Smith @*/ 463*bbf0e169SBarry Smith int MatGetInfo(Mat mat,MatInfoType flag,MatInfo *info) 464*bbf0e169SBarry Smith { 465*bbf0e169SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 466*bbf0e169SBarry Smith if (!mat->ops.getinfo) SETERRQ(PETSC_ERR_SUP,"MatGetInfo"); 467*bbf0e169SBarry Smith return (*mat->ops.getinfo)(mat,flag,info); 468*bbf0e169SBarry Smith } 469*bbf0e169SBarry Smith /* ----------------------------------------------------------*/ 470*bbf0e169SBarry Smith /*@ 471*bbf0e169SBarry Smith MatILUDTFactor - Performs a drop tolerance ILU factorization. 472*bbf0e169SBarry Smith 473*bbf0e169SBarry Smith Input Parameters: 474*bbf0e169SBarry Smith . mat - the matrix 475*bbf0e169SBarry Smith . dt - the drop tolerance 476*bbf0e169SBarry Smith . maxnz - the maximum number of nonzeros per row allowed? 477*bbf0e169SBarry Smith . row - row permutation 478*bbf0e169SBarry Smith . col - column permutation 479*bbf0e169SBarry Smith 480*bbf0e169SBarry Smith Output Parameters: 481*bbf0e169SBarry Smith . fact - the factored matrix 482*bbf0e169SBarry Smith 483*bbf0e169SBarry Smith .keywords: matrix, factor, LU, in-place 484*bbf0e169SBarry Smith 485*bbf0e169SBarry Smith .seealso: MatLUFactorSymbolic(), MatLUFactorNumeric(), MatCholeskyFactor() 486*bbf0e169SBarry Smith @*/ 487*bbf0e169SBarry Smith int MatILUDTFactor(Mat mat,double dt,int maxnz,IS row,IS col,Mat *fact) 488*bbf0e169SBarry Smith { 489*bbf0e169SBarry Smith int ierr; 490*bbf0e169SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 491*bbf0e169SBarry Smith if (!mat->ops.iludtfactor) SETERRQ(PETSC_ERR_SUP,"MatILUDTFactor"); 492*bbf0e169SBarry Smith if (!mat->assembled) SETERRQ(1,"MatILUDTFactor:Not for unassembled matrix"); 493*bbf0e169SBarry Smith if (mat->factor) SETERRQ(1,"MatILUDTFactor:Not for factored matrix"); 494*bbf0e169SBarry Smith 495*bbf0e169SBarry Smith PLogEventBegin(MAT_ILUFactor,mat,row,col,0); 496*bbf0e169SBarry Smith ierr = (*mat->ops.iludtfactor)(mat,dt,maxnz,row,col,fact); CHKERRQ(ierr); 497*bbf0e169SBarry Smith PLogEventEnd(MAT_ILUFactor,mat,row,col,0); 498*bbf0e169SBarry Smith 499*bbf0e169SBarry Smith return 0; 500*bbf0e169SBarry Smith } 501*bbf0e169SBarry Smith 502*bbf0e169SBarry Smith /*@ 503*bbf0e169SBarry Smith MatLUFactor - Performs in-place LU factorization of matrix. 504*bbf0e169SBarry Smith 505*bbf0e169SBarry Smith Input Parameters: 506*bbf0e169SBarry Smith . mat - the matrix 507*bbf0e169SBarry Smith . row - row permutation 508*bbf0e169SBarry Smith . col - column permutation 509*bbf0e169SBarry Smith . f - expected fill as ratio of original fill. 510*bbf0e169SBarry Smith 511*bbf0e169SBarry Smith .keywords: matrix, factor, LU, in-place 512*bbf0e169SBarry Smith 513*bbf0e169SBarry Smith .seealso: MatLUFactorSymbolic(), MatLUFactorNumeric(), MatCholeskyFactor() 514*bbf0e169SBarry Smith @*/ 515*bbf0e169SBarry Smith int MatLUFactor(Mat mat,IS row,IS col,double f) 516*bbf0e169SBarry Smith { 517*bbf0e169SBarry Smith int ierr; 518*bbf0e169SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 519*bbf0e169SBarry Smith if (mat->M != mat->N) SETERRQ(1,"MatLUFactor:matrix must be square"); 520*bbf0e169SBarry Smith if (!mat->ops.lufactor) SETERRQ(PETSC_ERR_SUP,"MatLUFactor"); 521*bbf0e169SBarry Smith if (!mat->assembled) SETERRQ(1,"MatLUFactor:Not for unassembled matrix"); 522*bbf0e169SBarry Smith if (mat->factor) SETERRQ(1,"MatLUFactor:Not for factored matrix"); 523*bbf0e169SBarry Smith 524*bbf0e169SBarry Smith PLogEventBegin(MAT_LUFactor,mat,row,col,0); 525*bbf0e169SBarry Smith ierr = (*mat->ops.lufactor)(mat,row,col,f); CHKERRQ(ierr); 526*bbf0e169SBarry Smith PLogEventEnd(MAT_LUFactor,mat,row,col,0); 527*bbf0e169SBarry Smith return 0; 528*bbf0e169SBarry Smith } 529*bbf0e169SBarry Smith /*@ 530*bbf0e169SBarry Smith MatILUFactor - Performs in-place ILU factorization of matrix. 531*bbf0e169SBarry Smith 532*bbf0e169SBarry Smith Input Parameters: 533*bbf0e169SBarry Smith . mat - the matrix 534*bbf0e169SBarry Smith . row - row permutation 535*bbf0e169SBarry Smith . col - column permutation 536*bbf0e169SBarry Smith . f - expected fill as ratio of original fill. 537*bbf0e169SBarry Smith . level - number of levels of fill. 538*bbf0e169SBarry Smith 539*bbf0e169SBarry Smith Note: probably really only in-place when level is zero. 540*bbf0e169SBarry Smith .keywords: matrix, factor, ILU, in-place 541*bbf0e169SBarry Smith 542*bbf0e169SBarry Smith .seealso: MatILUFactorSymbolic(), MatLUFactorNumeric(), MatCholeskyFactor() 543*bbf0e169SBarry Smith @*/ 544*bbf0e169SBarry Smith int MatILUFactor(Mat mat,IS row,IS col,double f,int level) 545*bbf0e169SBarry Smith { 546*bbf0e169SBarry Smith int ierr; 547*bbf0e169SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 548*bbf0e169SBarry Smith if (mat->M != mat->N) SETERRQ(1,"MatILUFactor:matrix must be square"); 549*bbf0e169SBarry Smith if (!mat->ops.ilufactor) SETERRQ(PETSC_ERR_SUP,"MatILUFactor"); 550*bbf0e169SBarry Smith if (!mat->assembled) SETERRQ(1,"MatILUFactor:Not for unassembled matrix"); 551*bbf0e169SBarry Smith if (mat->factor) SETERRQ(1,"MatILUFactor:Not for factored matrix"); 552*bbf0e169SBarry Smith 553*bbf0e169SBarry Smith PLogEventBegin(MAT_ILUFactor,mat,row,col,0); 554*bbf0e169SBarry Smith ierr = (*mat->ops.ilufactor)(mat,row,col,f,level); CHKERRQ(ierr); 555*bbf0e169SBarry Smith PLogEventEnd(MAT_ILUFactor,mat,row,col,0); 556*bbf0e169SBarry Smith return 0; 557*bbf0e169SBarry Smith } 558*bbf0e169SBarry Smith 559*bbf0e169SBarry Smith /*@ 560*bbf0e169SBarry Smith MatLUFactorSymbolic - Performs symbolic LU factorization of matrix. 561*bbf0e169SBarry Smith Call this routine before calling MatLUFactorNumeric(). 562*bbf0e169SBarry Smith 563*bbf0e169SBarry Smith Input Parameters: 564*bbf0e169SBarry Smith . mat - the matrix 565*bbf0e169SBarry Smith . row, col - row and column permutations 566*bbf0e169SBarry Smith . f - expected fill as ratio of the original number of nonzeros, 567*bbf0e169SBarry Smith for example 3.0; choosing this parameter well can result in 568*bbf0e169SBarry Smith more efficient use of time and space. 569*bbf0e169SBarry Smith 570*bbf0e169SBarry Smith Output Parameter: 571*bbf0e169SBarry Smith . fact - new matrix that has been symbolically factored 572*bbf0e169SBarry Smith 573*bbf0e169SBarry Smith Options Database Key: 574*bbf0e169SBarry Smith $ -mat_lu_fill <f>, where f is the fill ratio 575*bbf0e169SBarry Smith 576*bbf0e169SBarry Smith Notes: 577*bbf0e169SBarry Smith See the file $(PETSC_DIR)/Performace for additional information about 578*bbf0e169SBarry Smith choosing the fill factor for better efficiency. 579*bbf0e169SBarry Smith 580*bbf0e169SBarry Smith .keywords: matrix, factor, LU, symbolic, fill 581*bbf0e169SBarry Smith 582*bbf0e169SBarry Smith .seealso: MatLUFactor(), MatLUFactorNumeric(), MatCholeskyFactor() 583*bbf0e169SBarry Smith @*/ 584*bbf0e169SBarry Smith int MatLUFactorSymbolic(Mat mat,IS row,IS col,double f,Mat *fact) 585*bbf0e169SBarry Smith { 586*bbf0e169SBarry Smith int ierr,flg; 587*bbf0e169SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 588*bbf0e169SBarry Smith if (mat->M != mat->N) SETERRQ(1,"MatLUFactorSymbolic:matrix must be square"); 589*bbf0e169SBarry Smith if (!fact) SETERRQ(1,"MatLUFactorSymbolic:Missing factor matrix argument"); 590*bbf0e169SBarry Smith if (!mat->ops.lufactorsymbolic) SETERRQ(PETSC_ERR_SUP,"MatLUFactorSymbolic"); 591*bbf0e169SBarry Smith if (!mat->assembled) SETERRQ(1,"MatLUFactorSymbolic:Not for unassembled matrix"); 592*bbf0e169SBarry Smith if (mat->factor) SETERRQ(1,"MatLUFactorSymbolic:Not for factored matrix"); 593*bbf0e169SBarry Smith 594*bbf0e169SBarry Smith ierr = OptionsGetDouble(PETSC_NULL,"-mat_lu_fill",&f,&flg); CHKERRQ(ierr); 595*bbf0e169SBarry Smith PLogEventBegin(MAT_LUFactorSymbolic,mat,row,col,0); 596*bbf0e169SBarry Smith ierr = (*mat->ops.lufactorsymbolic)(mat,row,col,f,fact); CHKERRQ(ierr); 597*bbf0e169SBarry Smith PLogEventEnd(MAT_LUFactorSymbolic,mat,row,col,0); 598*bbf0e169SBarry Smith return 0; 599*bbf0e169SBarry Smith } 600*bbf0e169SBarry Smith /*@ 601*bbf0e169SBarry Smith MatLUFactorNumeric - Performs numeric LU factorization of a matrix. 602*bbf0e169SBarry Smith Call this routine after first calling MatLUFactorSymbolic(). 603*bbf0e169SBarry Smith 604*bbf0e169SBarry Smith Input Parameters: 605*bbf0e169SBarry Smith . mat - the matrix 606*bbf0e169SBarry Smith . row, col - row and column permutations 607*bbf0e169SBarry Smith 608*bbf0e169SBarry Smith Output Parameters: 609*bbf0e169SBarry Smith . fact - symbolically factored matrix that must have been generated 610*bbf0e169SBarry Smith by MatLUFactorSymbolic() 611*bbf0e169SBarry Smith 612*bbf0e169SBarry Smith Notes: 613*bbf0e169SBarry Smith See MatLUFactor() for in-place factorization. See 614*bbf0e169SBarry Smith MatCholeskyFactorNumeric() for the symmetric, positive definite case. 615*bbf0e169SBarry Smith 616*bbf0e169SBarry Smith .keywords: matrix, factor, LU, numeric 617*bbf0e169SBarry Smith 618*bbf0e169SBarry Smith .seealso: MatLUFactorSymbolic(), MatLUFactor(), MatCholeskyFactor() 619*bbf0e169SBarry Smith @*/ 620*bbf0e169SBarry Smith int MatLUFactorNumeric(Mat mat,Mat *fact) 621*bbf0e169SBarry Smith { 622*bbf0e169SBarry Smith int ierr,flg; 623*bbf0e169SBarry Smith 624*bbf0e169SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 625*bbf0e169SBarry Smith if (!fact) SETERRQ(1,"MatLUFactorNumeric:Missing factor matrix argument"); 626*bbf0e169SBarry Smith if (!mat->ops.lufactornumeric) SETERRQ(PETSC_ERR_SUP,"MatLUFactorNumeric"); 627*bbf0e169SBarry Smith if (!mat->assembled) SETERRQ(1,"MatLUFactorNumeric:Not for unassembled matrix"); 628*bbf0e169SBarry Smith if (mat->M != (*fact)->M || mat->N != (*fact)->N) 629*bbf0e169SBarry Smith SETERRQ(PETSC_ERR_SIZ,"MatLUFactorNumeric:Mat mat,Mat *fact: global dim"); 630*bbf0e169SBarry Smith 631*bbf0e169SBarry Smith PLogEventBegin(MAT_LUFactorNumeric,mat,*fact,0,0); 632*bbf0e169SBarry Smith ierr = (*mat->ops.lufactornumeric)(mat,fact); CHKERRQ(ierr); 633*bbf0e169SBarry Smith PLogEventEnd(MAT_LUFactorNumeric,mat,*fact,0,0); 634*bbf0e169SBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_view_draw",&flg); CHKERRQ(ierr); 635*bbf0e169SBarry Smith if (flg) { 636*bbf0e169SBarry Smith Viewer viewer; 637*bbf0e169SBarry Smith ierr = ViewerDrawOpenX((*fact)->comm,0,0,0,0,300,300,&viewer);CHKERRQ(ierr); 638*bbf0e169SBarry Smith ierr = MatView(*fact,viewer); CHKERRQ(ierr); 639*bbf0e169SBarry Smith ierr = ViewerFlush(viewer); CHKERRQ(ierr); 640*bbf0e169SBarry Smith ierr = ViewerDestroy(viewer); CHKERRQ(ierr); 641*bbf0e169SBarry Smith } 642*bbf0e169SBarry Smith return 0; 643*bbf0e169SBarry Smith } 644*bbf0e169SBarry Smith /*@ 645*bbf0e169SBarry Smith MatCholeskyFactor - Performs in-place Cholesky factorization of a 646*bbf0e169SBarry Smith symmetric matrix. 647*bbf0e169SBarry Smith 648*bbf0e169SBarry Smith Input Parameters: 649*bbf0e169SBarry Smith . mat - the matrix 650*bbf0e169SBarry Smith . perm - row and column permutations 651*bbf0e169SBarry Smith . f - expected fill as ratio of original fill 652*bbf0e169SBarry Smith 653*bbf0e169SBarry Smith Notes: 654*bbf0e169SBarry Smith See MatLUFactor() for the nonsymmetric case. See also 655*bbf0e169SBarry Smith MatCholeskyFactorSymbolic(), and MatCholeskyFactorNumeric(). 656*bbf0e169SBarry Smith 657*bbf0e169SBarry Smith .keywords: matrix, factor, in-place, Cholesky 658*bbf0e169SBarry Smith 659*bbf0e169SBarry Smith .seealso: MatLUFactor(), MatCholeskyFactorSymbolic(), MatCholeskyFactorNumeric() 660*bbf0e169SBarry Smith @*/ 661*bbf0e169SBarry Smith int MatCholeskyFactor(Mat mat,IS perm,double f) 662*bbf0e169SBarry Smith { 663*bbf0e169SBarry Smith int ierr; 664*bbf0e169SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 665*bbf0e169SBarry Smith if (mat->M != mat->N) SETERRQ(1,"MatCholeskyFactor:matrix must be square"); 666*bbf0e169SBarry Smith if (!mat->ops.choleskyfactor) SETERRQ(PETSC_ERR_SUP,"MatCholeskyFactor"); 667*bbf0e169SBarry Smith if (!mat->assembled) SETERRQ(1,"MatCholeskyFactor:Not for unassembled matrix"); 668*bbf0e169SBarry Smith if (mat->factor) SETERRQ(1,"MatCholeskyFactor:Not for factored matrix"); 669*bbf0e169SBarry Smith 670*bbf0e169SBarry Smith PLogEventBegin(MAT_CholeskyFactor,mat,perm,0,0); 671*bbf0e169SBarry Smith ierr = (*mat->ops.choleskyfactor)(mat,perm,f); CHKERRQ(ierr); 672*bbf0e169SBarry Smith PLogEventEnd(MAT_CholeskyFactor,mat,perm,0,0); 673*bbf0e169SBarry Smith return 0; 674*bbf0e169SBarry Smith } 675*bbf0e169SBarry Smith /*@ 676*bbf0e169SBarry Smith MatCholeskyFactorSymbolic - Performs symbolic Cholesky factorization 677*bbf0e169SBarry Smith of a symmetric matrix. 678*bbf0e169SBarry Smith 679*bbf0e169SBarry Smith Input Parameters: 680*bbf0e169SBarry Smith . mat - the matrix 681*bbf0e169SBarry Smith . perm - row and column permutations 682*bbf0e169SBarry Smith . f - expected fill as ratio of original 683*bbf0e169SBarry Smith 684*bbf0e169SBarry Smith Output Parameter: 685*bbf0e169SBarry Smith . fact - the factored matrix 686*bbf0e169SBarry Smith 687*bbf0e169SBarry Smith Notes: 688*bbf0e169SBarry Smith See MatLUFactorSymbolic() for the nonsymmetric case. See also 689*bbf0e169SBarry Smith MatCholeskyFactor() and MatCholeskyFactorNumeric(). 690*bbf0e169SBarry Smith 691*bbf0e169SBarry Smith .keywords: matrix, factor, factorization, symbolic, Cholesky 692*bbf0e169SBarry Smith 693*bbf0e169SBarry Smith .seealso: MatLUFactorSymbolic(), MatCholeskyFactor(), MatCholeskyFactorNumeric() 694*bbf0e169SBarry Smith @*/ 695*bbf0e169SBarry Smith int MatCholeskyFactorSymbolic(Mat mat,IS perm,double f,Mat *fact) 696*bbf0e169SBarry Smith { 697*bbf0e169SBarry Smith int ierr; 698*bbf0e169SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 699*bbf0e169SBarry Smith if (mat->M != mat->N) SETERRQ(1,"MatCholeskyFactorSymbolic:matrix must be square"); 700*bbf0e169SBarry Smith if (!fact) SETERRQ(1,"MatCholeskyFactorSymbolic:Missing factor matrix argument"); 701*bbf0e169SBarry Smith if (!mat->ops.choleskyfactorsymbolic)SETERRQ(PETSC_ERR_SUP,"MatCholeskyFactorSymbolic"); 702*bbf0e169SBarry Smith if (!mat->assembled) SETERRQ(1,"MatCholeskyFactorSymbolic:Not for unassembled matrix"); 703*bbf0e169SBarry Smith if (mat->factor) SETERRQ(1,"MatCholeskyFactorSymbolic:Not for factored matrix"); 704*bbf0e169SBarry Smith 705*bbf0e169SBarry Smith PLogEventBegin(MAT_CholeskyFactorSymbolic,mat,perm,0,0); 706*bbf0e169SBarry Smith ierr = (*mat->ops.choleskyfactorsymbolic)(mat,perm,f,fact); CHKERRQ(ierr); 707*bbf0e169SBarry Smith PLogEventEnd(MAT_CholeskyFactorSymbolic,mat,perm,0,0); 708*bbf0e169SBarry Smith return 0; 709*bbf0e169SBarry Smith } 710*bbf0e169SBarry Smith /*@ 711*bbf0e169SBarry Smith MatCholeskyFactorNumeric - Performs numeric Cholesky factorization 712*bbf0e169SBarry Smith of a symmetric matrix. Call this routine after first calling 713*bbf0e169SBarry Smith MatCholeskyFactorSymbolic(). 714*bbf0e169SBarry Smith 715*bbf0e169SBarry Smith Input Parameter: 716*bbf0e169SBarry Smith . mat - the initial matrix 717*bbf0e169SBarry Smith 718*bbf0e169SBarry Smith Output Parameter: 719*bbf0e169SBarry Smith . fact - the factored matrix 720*bbf0e169SBarry Smith 721*bbf0e169SBarry Smith .keywords: matrix, factor, numeric, Cholesky 722*bbf0e169SBarry Smith 723*bbf0e169SBarry Smith .seealso: MatCholeskyFactorSymbolic(), MatCholeskyFactor(), MatLUFactorNumeric() 724*bbf0e169SBarry Smith @*/ 725*bbf0e169SBarry Smith int MatCholeskyFactorNumeric(Mat mat,Mat *fact) 726*bbf0e169SBarry Smith { 727*bbf0e169SBarry Smith int ierr; 728*bbf0e169SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 729*bbf0e169SBarry Smith if (!fact) SETERRQ(1,"MatCholeskyFactorNumeric:Missing factor matrix argument"); 730*bbf0e169SBarry Smith if (!mat->ops.choleskyfactornumeric) SETERRQ(PETSC_ERR_SUP,"MatCholeskyFactorNumeric"); 731*bbf0e169SBarry Smith if (!mat->assembled) SETERRQ(1,"MatCholeskyFactorNumeric:Not for unassembled matrix"); 732*bbf0e169SBarry Smith if (mat->M != (*fact)->M || mat->N != (*fact)->N) 733*bbf0e169SBarry Smith SETERRQ(PETSC_ERR_SIZ,"MatCholeskyFactorNumeric:Mat mat,Mat *fact: global dim"); 734*bbf0e169SBarry Smith 735*bbf0e169SBarry Smith PLogEventBegin(MAT_CholeskyFactorNumeric,mat,*fact,0,0); 736*bbf0e169SBarry Smith ierr = (*mat->ops.choleskyfactornumeric)(mat,fact); CHKERRQ(ierr); 737*bbf0e169SBarry Smith PLogEventEnd(MAT_CholeskyFactorNumeric,mat,*fact,0,0); 738*bbf0e169SBarry Smith return 0; 739*bbf0e169SBarry Smith } 740*bbf0e169SBarry Smith /* ----------------------------------------------------------------*/ 741*bbf0e169SBarry Smith /*@ 742*bbf0e169SBarry Smith MatSolve - Solves A x = b, given a factored matrix. 743*bbf0e169SBarry Smith 744*bbf0e169SBarry Smith Input Parameters: 745*bbf0e169SBarry Smith . mat - the factored matrix 746*bbf0e169SBarry Smith . b - the right-hand-side vector 747*bbf0e169SBarry Smith 748*bbf0e169SBarry Smith Output Parameter: 749*bbf0e169SBarry Smith . x - the result vector 750*bbf0e169SBarry Smith 751*bbf0e169SBarry Smith .keywords: matrix, linear system, solve, LU, Cholesky, triangular solve 752*bbf0e169SBarry Smith 753*bbf0e169SBarry Smith .seealso: MatSolveAdd(), MatSolveTrans(), MatSolveTransAdd() 754*bbf0e169SBarry Smith @*/ 755*bbf0e169SBarry Smith int MatSolve(Mat mat,Vec b,Vec x) 756*bbf0e169SBarry Smith { 757*bbf0e169SBarry Smith int ierr; 758*bbf0e169SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 759*bbf0e169SBarry Smith PetscValidHeaderSpecific(b,VEC_COOKIE); PetscValidHeaderSpecific(x,VEC_COOKIE); 760*bbf0e169SBarry Smith if (x == b) SETERRQ(1,"MatSolve:x and y must be different vectors"); 761*bbf0e169SBarry Smith if (!mat->factor) SETERRQ(1,"MatSolve:Unfactored matrix"); 762*bbf0e169SBarry Smith if (mat->N != x->N) SETERRQ(PETSC_ERR_SIZ,"MatSolve:Mat mat,Vec x: global dim"); 763*bbf0e169SBarry Smith if (mat->M != b->N) SETERRQ(PETSC_ERR_SIZ,"MatSolve:Mat mat,Vec b: global dim"); 764*bbf0e169SBarry Smith if (mat->m != b->n) SETERRQ(PETSC_ERR_SIZ,"MatSolve:Mat mat,Vec b: local dim"); 765*bbf0e169SBarry Smith 766*bbf0e169SBarry Smith if (!mat->ops.solve) SETERRQ(PETSC_ERR_SUP,"MatSolve"); 767*bbf0e169SBarry Smith PLogEventBegin(MAT_Solve,mat,b,x,0); 768*bbf0e169SBarry Smith ierr = (*mat->ops.solve)(mat,b,x); CHKERRQ(ierr); 769*bbf0e169SBarry Smith PLogEventEnd(MAT_Solve,mat,b,x,0); 770*bbf0e169SBarry Smith return 0; 771*bbf0e169SBarry Smith } 772*bbf0e169SBarry Smith 773*bbf0e169SBarry Smith /* @ 774*bbf0e169SBarry Smith MatForwardSolve - Solves L x = b, given a factored matrix, A = LU. 775*bbf0e169SBarry Smith 776*bbf0e169SBarry Smith Input Parameters: 777*bbf0e169SBarry Smith . mat - the factored matrix 778*bbf0e169SBarry Smith . b - the right-hand-side vector 779*bbf0e169SBarry Smith 780*bbf0e169SBarry Smith Output Parameter: 781*bbf0e169SBarry Smith . x - the result vector 782*bbf0e169SBarry Smith 783*bbf0e169SBarry Smith Notes: 784*bbf0e169SBarry Smith MatSolve() should be used for most applications, as it performs 785*bbf0e169SBarry Smith a forward solve followed by a backward solve. 786*bbf0e169SBarry Smith 787*bbf0e169SBarry Smith .keywords: matrix, forward, LU, Cholesky, triangular solve 788*bbf0e169SBarry Smith 789*bbf0e169SBarry Smith .seealso: MatSolve(), MatBackwardSolve() 790*bbf0e169SBarry Smith @ */ 791*bbf0e169SBarry Smith int MatForwardSolve(Mat mat,Vec b,Vec x) 792*bbf0e169SBarry Smith { 793*bbf0e169SBarry Smith int ierr; 794*bbf0e169SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 795*bbf0e169SBarry Smith PetscValidHeaderSpecific(b,VEC_COOKIE); PetscValidHeaderSpecific(x,VEC_COOKIE); 796*bbf0e169SBarry Smith if (x == b) SETERRQ(1,"MatForwardSolve:x and y must be different vectors"); 797*bbf0e169SBarry Smith if (!mat->factor) SETERRQ(1,"MatForwardSolve:Unfactored matrix"); 798*bbf0e169SBarry Smith if (!mat->ops.forwardsolve) SETERRQ(PETSC_ERR_SUP,"MatForwardSolve"); 799*bbf0e169SBarry Smith if (mat->N != x->N) SETERRQ(PETSC_ERR_SIZ,"MatForwardSolve:Mat mat,Vec x: global dim"); 800*bbf0e169SBarry Smith if (mat->M != b->N) SETERRQ(PETSC_ERR_SIZ,"MatForwardSolve:Mat mat,Vec b: global dim"); 801*bbf0e169SBarry Smith if (mat->m != b->n) SETERRQ(PETSC_ERR_SIZ,"MatForwardSolve:Mat mat,Vec b: local dim"); 802*bbf0e169SBarry Smith 803*bbf0e169SBarry Smith PLogEventBegin(MAT_ForwardSolve,mat,b,x,0); 804*bbf0e169SBarry Smith ierr = (*mat->ops.forwardsolve)(mat,b,x); CHKERRQ(ierr); 805*bbf0e169SBarry Smith PLogEventEnd(MAT_ForwardSolve,mat,b,x,0); 806*bbf0e169SBarry Smith return 0; 807*bbf0e169SBarry Smith } 808*bbf0e169SBarry Smith 809*bbf0e169SBarry Smith /* @ 810*bbf0e169SBarry Smith MatBackwardSolve - Solves U x = b, given a factored matrix, A = LU. 811*bbf0e169SBarry Smith 812*bbf0e169SBarry Smith Input Parameters: 813*bbf0e169SBarry Smith . mat - the factored matrix 814*bbf0e169SBarry Smith . b - the right-hand-side vector 815*bbf0e169SBarry Smith 816*bbf0e169SBarry Smith Output Parameter: 817*bbf0e169SBarry Smith . x - the result vector 818*bbf0e169SBarry Smith 819*bbf0e169SBarry Smith Notes: 820*bbf0e169SBarry Smith MatSolve() should be used for most applications, as it performs 821*bbf0e169SBarry Smith a forward solve followed by a backward solve. 822*bbf0e169SBarry Smith 823*bbf0e169SBarry Smith .keywords: matrix, backward, LU, Cholesky, triangular solve 824*bbf0e169SBarry Smith 825*bbf0e169SBarry Smith .seealso: MatSolve(), MatForwardSolve() 826*bbf0e169SBarry Smith @ */ 827*bbf0e169SBarry Smith int MatBackwardSolve(Mat mat,Vec b,Vec x) 828*bbf0e169SBarry Smith { 829*bbf0e169SBarry Smith int ierr; 830*bbf0e169SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 831*bbf0e169SBarry Smith PetscValidHeaderSpecific(b,VEC_COOKIE); PetscValidHeaderSpecific(x,VEC_COOKIE); 832*bbf0e169SBarry Smith if (x == b) SETERRQ(1,"MatBackwardSolve:x and b must be different vectors"); 833*bbf0e169SBarry Smith if (!mat->factor) SETERRQ(1,"MatBackwardSolve:Unfactored matrix"); 834*bbf0e169SBarry Smith if (!mat->ops.backwardsolve) SETERRQ(PETSC_ERR_SUP,"MatBackwardSolve"); 835*bbf0e169SBarry Smith if (mat->N != x->N) SETERRQ(PETSC_ERR_SIZ,"MatBackwardSolve:Mat mat,Vec x: global dim"); 836*bbf0e169SBarry Smith if (mat->M != b->N) SETERRQ(PETSC_ERR_SIZ,"MatBackwardSolve:Mat mat,Vec b: global dim"); 837*bbf0e169SBarry Smith if (mat->m != b->n) SETERRQ(PETSC_ERR_SIZ,"MatBackwardSolve:Mat mat,Vec b: local dim"); 838*bbf0e169SBarry Smith 839*bbf0e169SBarry Smith PLogEventBegin(MAT_BackwardSolve,mat,b,x,0); 840*bbf0e169SBarry Smith ierr = (*mat->ops.backwardsolve)(mat,b,x); CHKERRQ(ierr); 841*bbf0e169SBarry Smith PLogEventEnd(MAT_BackwardSolve,mat,b,x,0); 842*bbf0e169SBarry Smith return 0; 843*bbf0e169SBarry Smith } 844*bbf0e169SBarry Smith 845*bbf0e169SBarry Smith /*@ 846*bbf0e169SBarry Smith MatSolveAdd - Computes x = y + inv(A)*b, given a factored matrix. 847*bbf0e169SBarry Smith 848*bbf0e169SBarry Smith Input Parameters: 849*bbf0e169SBarry Smith . mat - the factored matrix 850*bbf0e169SBarry Smith . b - the right-hand-side vector 851*bbf0e169SBarry Smith . y - the vector to be added to 852*bbf0e169SBarry Smith 853*bbf0e169SBarry Smith Output Parameter: 854*bbf0e169SBarry Smith . x - the result vector 855*bbf0e169SBarry Smith 856*bbf0e169SBarry Smith .keywords: matrix, linear system, solve, LU, Cholesky, add 857*bbf0e169SBarry Smith 858*bbf0e169SBarry Smith .seealso: MatSolve(), MatSolveTrans(), MatSolveTransAdd() 859*bbf0e169SBarry Smith @*/ 860*bbf0e169SBarry Smith int MatSolveAdd(Mat mat,Vec b,Vec y,Vec x) 861*bbf0e169SBarry Smith { 862*bbf0e169SBarry Smith Scalar one = 1.0; 863*bbf0e169SBarry Smith Vec tmp; 864*bbf0e169SBarry Smith int ierr; 865*bbf0e169SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE);PetscValidHeaderSpecific(y,VEC_COOKIE); 866*bbf0e169SBarry Smith PetscValidHeaderSpecific(b,VEC_COOKIE); PetscValidHeaderSpecific(x,VEC_COOKIE); 867*bbf0e169SBarry Smith if (x == b) SETERRQ(1,"MatSolveAdd:x and b must be different vectors"); 868*bbf0e169SBarry Smith if (!mat->factor) SETERRQ(1,"MatSolveAdd:Unfactored matrix"); 869*bbf0e169SBarry Smith if (mat->N != x->N) SETERRQ(PETSC_ERR_SIZ,"MatSolveAdd:Mat mat,Vec x: global dim"); 870*bbf0e169SBarry Smith if (mat->M != b->N) SETERRQ(PETSC_ERR_SIZ,"MatSolveAdd:Mat mat,Vec b: global dim"); 871*bbf0e169SBarry Smith if (mat->M != y->N) SETERRQ(PETSC_ERR_SIZ,"MatSolveAdd:Mat mat,Vec y: global dim"); 872*bbf0e169SBarry Smith if (mat->m != b->n) SETERRQ(PETSC_ERR_SIZ,"MatSolveAdd:Mat mat,Vec b: local dim"); 873*bbf0e169SBarry Smith if (x->n != y->n) SETERRQ(PETSC_ERR_SIZ,"MatSolveAdd:Vec x,Vec y: local dim"); 874*bbf0e169SBarry Smith 875*bbf0e169SBarry Smith PLogEventBegin(MAT_SolveAdd,mat,b,x,y); 876*bbf0e169SBarry Smith if (mat->ops.solveadd) { 877*bbf0e169SBarry Smith ierr = (*mat->ops.solveadd)(mat,b,y,x); CHKERRQ(ierr); 878*bbf0e169SBarry Smith } 879*bbf0e169SBarry Smith else { 880*bbf0e169SBarry Smith /* do the solve then the add manually */ 881*bbf0e169SBarry Smith if (x != y) { 882*bbf0e169SBarry Smith ierr = MatSolve(mat,b,x); CHKERRQ(ierr); 883*bbf0e169SBarry Smith ierr = VecAXPY(&one,y,x); CHKERRQ(ierr); 884*bbf0e169SBarry Smith } 885*bbf0e169SBarry Smith else { 886*bbf0e169SBarry Smith ierr = VecDuplicate(x,&tmp); CHKERRQ(ierr); 887*bbf0e169SBarry Smith PLogObjectParent(mat,tmp); 888*bbf0e169SBarry Smith ierr = VecCopy(x,tmp); CHKERRQ(ierr); 889*bbf0e169SBarry Smith ierr = MatSolve(mat,b,x); CHKERRQ(ierr); 890*bbf0e169SBarry Smith ierr = VecAXPY(&one,tmp,x); CHKERRQ(ierr); 891*bbf0e169SBarry Smith ierr = VecDestroy(tmp); CHKERRQ(ierr); 892*bbf0e169SBarry Smith } 893*bbf0e169SBarry Smith } 894*bbf0e169SBarry Smith PLogEventEnd(MAT_SolveAdd,mat,b,x,y); 895*bbf0e169SBarry Smith return 0; 896*bbf0e169SBarry Smith } 897*bbf0e169SBarry Smith /*@ 898*bbf0e169SBarry Smith MatSolveTrans - Solves A' x = b, given a factored matrix. 899*bbf0e169SBarry Smith 900*bbf0e169SBarry Smith Input Parameters: 901*bbf0e169SBarry Smith . mat - the factored matrix 902*bbf0e169SBarry Smith . b - the right-hand-side vector 903*bbf0e169SBarry Smith 904*bbf0e169SBarry Smith Output Parameter: 905*bbf0e169SBarry Smith . x - the result vector 906*bbf0e169SBarry Smith 907*bbf0e169SBarry Smith .keywords: matrix, linear system, solve, LU, Cholesky, transpose 908*bbf0e169SBarry Smith 909*bbf0e169SBarry Smith .seealso: MatSolve(), MatSolveAdd(), MatSolveTransAdd() 910*bbf0e169SBarry Smith @*/ 911*bbf0e169SBarry Smith int MatSolveTrans(Mat mat,Vec b,Vec x) 912*bbf0e169SBarry Smith { 913*bbf0e169SBarry Smith int ierr; 914*bbf0e169SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 915*bbf0e169SBarry Smith PetscValidHeaderSpecific(b,VEC_COOKIE); PetscValidHeaderSpecific(x,VEC_COOKIE); 916*bbf0e169SBarry Smith if (!mat->factor) SETERRQ(1,"MatSolveTrans:Unfactored matrix"); 917*bbf0e169SBarry Smith if (x == b) SETERRQ(1,"MatSolveTrans:x and b must be different vectors"); 918*bbf0e169SBarry Smith if (!mat->ops.solvetrans) SETERRQ(PETSC_ERR_SUP,"MatSolveTrans"); 919*bbf0e169SBarry Smith if (mat->M != x->N) SETERRQ(PETSC_ERR_SIZ,"MatSolveTrans:Mat mat,Vec x: global dim"); 920*bbf0e169SBarry Smith if (mat->N != b->N) SETERRQ(PETSC_ERR_SIZ,"MatSolveTrans:Mat mat,Vec b: global dim"); 921*bbf0e169SBarry Smith 922*bbf0e169SBarry Smith PLogEventBegin(MAT_SolveTrans,mat,b,x,0); 923*bbf0e169SBarry Smith ierr = (*mat->ops.solvetrans)(mat,b,x); CHKERRQ(ierr); 924*bbf0e169SBarry Smith PLogEventEnd(MAT_SolveTrans,mat,b,x,0); 925*bbf0e169SBarry Smith return 0; 926*bbf0e169SBarry Smith } 927*bbf0e169SBarry Smith /*@ 928*bbf0e169SBarry Smith MatSolveTransAdd - Computes x = y + inv(trans(A)) b, given a 929*bbf0e169SBarry Smith factored matrix. 930*bbf0e169SBarry Smith 931*bbf0e169SBarry Smith Input Parameters: 932*bbf0e169SBarry Smith . mat - the factored matrix 933*bbf0e169SBarry Smith . b - the right-hand-side vector 934*bbf0e169SBarry Smith . y - the vector to be added to 935*bbf0e169SBarry Smith 936*bbf0e169SBarry Smith Output Parameter: 937*bbf0e169SBarry Smith . x - the result vector 938*bbf0e169SBarry Smith 939*bbf0e169SBarry Smith .keywords: matrix, linear system, solve, LU, Cholesky, transpose, add 940*bbf0e169SBarry Smith 941*bbf0e169SBarry Smith .seealso: MatSolve(), MatSolveAdd(), MatSolveTrans() 942*bbf0e169SBarry Smith @*/ 943*bbf0e169SBarry Smith int MatSolveTransAdd(Mat mat,Vec b,Vec y,Vec x) 944*bbf0e169SBarry Smith { 945*bbf0e169SBarry Smith Scalar one = 1.0; 946*bbf0e169SBarry Smith int ierr; 947*bbf0e169SBarry Smith Vec tmp; 948*bbf0e169SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE);PetscValidHeaderSpecific(y,VEC_COOKIE); 949*bbf0e169SBarry Smith PetscValidHeaderSpecific(b,VEC_COOKIE); PetscValidHeaderSpecific(x,VEC_COOKIE); 950*bbf0e169SBarry Smith if (x == b) SETERRQ(1,"MatSolveTransAdd:x and b must be different vectors"); 951*bbf0e169SBarry Smith if (!mat->factor) SETERRQ(1,"MatSolveTransAdd:Unfactored matrix"); 952*bbf0e169SBarry Smith if (mat->M != x->N) SETERRQ(PETSC_ERR_SIZ,"MatSolveTransAdd:Mat mat,Vec x: global dim"); 953*bbf0e169SBarry Smith if (mat->N != b->N) SETERRQ(PETSC_ERR_SIZ,"MatSolveTransAdd:Mat mat,Vec b: global dim"); 954*bbf0e169SBarry Smith if (mat->N != y->N) SETERRQ(PETSC_ERR_SIZ,"MatSolveTransAdd:Mat mat,Vec y: global dim"); 955*bbf0e169SBarry Smith if (x->n != y->n) SETERRQ(PETSC_ERR_SIZ,"MatSolveTransAdd:Vec x,Vec y: local dim"); 956*bbf0e169SBarry Smith 957*bbf0e169SBarry Smith PLogEventBegin(MAT_SolveTransAdd,mat,b,x,y); 958*bbf0e169SBarry Smith if (mat->ops.solvetransadd) { 959*bbf0e169SBarry Smith ierr = (*mat->ops.solvetransadd)(mat,b,y,x); CHKERRQ(ierr); 960*bbf0e169SBarry Smith } 961*bbf0e169SBarry Smith else { 962*bbf0e169SBarry Smith /* do the solve then the add manually */ 963*bbf0e169SBarry Smith if (x != y) { 964*bbf0e169SBarry Smith ierr = MatSolveTrans(mat,b,x); CHKERRQ(ierr); 965*bbf0e169SBarry Smith ierr = VecAXPY(&one,y,x); CHKERRQ(ierr); 966*bbf0e169SBarry Smith } 967*bbf0e169SBarry Smith else { 968*bbf0e169SBarry Smith ierr = VecDuplicate(x,&tmp); CHKERRQ(ierr); 969*bbf0e169SBarry Smith PLogObjectParent(mat,tmp); 970*bbf0e169SBarry Smith ierr = VecCopy(x,tmp); CHKERRQ(ierr); 971*bbf0e169SBarry Smith ierr = MatSolveTrans(mat,b,x); CHKERRQ(ierr); 972*bbf0e169SBarry Smith ierr = VecAXPY(&one,tmp,x); CHKERRQ(ierr); 973*bbf0e169SBarry Smith ierr = VecDestroy(tmp); CHKERRQ(ierr); 974*bbf0e169SBarry Smith } 975*bbf0e169SBarry Smith } 976*bbf0e169SBarry Smith PLogEventEnd(MAT_SolveTransAdd,mat,b,x,y); 977*bbf0e169SBarry Smith return 0; 978*bbf0e169SBarry Smith } 979*bbf0e169SBarry Smith /* ----------------------------------------------------------------*/ 980*bbf0e169SBarry Smith 981*bbf0e169SBarry Smith /*@ 982*bbf0e169SBarry Smith MatRelax - Computes one relaxation sweep. 983*bbf0e169SBarry Smith 984*bbf0e169SBarry Smith Input Parameters: 985*bbf0e169SBarry Smith . mat - the matrix 986*bbf0e169SBarry Smith . b - the right hand side 987*bbf0e169SBarry Smith . omega - the relaxation factor 988*bbf0e169SBarry Smith . flag - flag indicating the type of SOR, one of 989*bbf0e169SBarry Smith $ SOR_FORWARD_SWEEP 990*bbf0e169SBarry Smith $ SOR_BACKWARD_SWEEP 991*bbf0e169SBarry Smith $ SOR_SYMMETRIC_SWEEP (SSOR method) 992*bbf0e169SBarry Smith $ SOR_LOCAL_FORWARD_SWEEP 993*bbf0e169SBarry Smith $ SOR_LOCAL_BACKWARD_SWEEP 994*bbf0e169SBarry Smith $ SOR_LOCAL_SYMMETRIC_SWEEP (local SSOR) 995*bbf0e169SBarry Smith $ SOR_APPLY_UPPER, SOR_APPLY_LOWER - applies 996*bbf0e169SBarry Smith $ upper/lower triangular part of matrix to 997*bbf0e169SBarry Smith $ vector (with omega) 998*bbf0e169SBarry Smith $ SOR_ZERO_INITIAL_GUESS - zero initial guess 999*bbf0e169SBarry Smith . shift - diagonal shift 1000*bbf0e169SBarry Smith . its - the number of iterations 1001*bbf0e169SBarry Smith 1002*bbf0e169SBarry Smith Output Parameters: 1003*bbf0e169SBarry Smith . x - the solution (can contain an initial guess) 1004*bbf0e169SBarry Smith 1005*bbf0e169SBarry Smith Notes: 1006*bbf0e169SBarry Smith SOR_LOCAL_FORWARD_SWEEP, SOR_LOCAL_BACKWARD_SWEEP, and 1007*bbf0e169SBarry Smith SOR_LOCAL_SYMMETRIC_SWEEP perform seperate independent smoothings 1008*bbf0e169SBarry Smith on each processor. 1009*bbf0e169SBarry Smith 1010*bbf0e169SBarry Smith Application programmers will not generally use MatRelax() directly, 1011*bbf0e169SBarry Smith but instead will employ the SLES/PC interface. 1012*bbf0e169SBarry Smith 1013*bbf0e169SBarry Smith Notes for Advanced Users: 1014*bbf0e169SBarry Smith The flags are implemented as bitwise inclusive or operations. 1015*bbf0e169SBarry Smith For example, use (SOR_ZERO_INITIAL_GUESS | SOR_SYMMETRIC_SWEEP) 1016*bbf0e169SBarry Smith to specify a zero initial guess for SSOR. 1017*bbf0e169SBarry Smith 1018*bbf0e169SBarry Smith .keywords: matrix, relax, relaxation, sweep 1019*bbf0e169SBarry Smith @*/ 1020*bbf0e169SBarry Smith int MatRelax(Mat mat,Vec b,double omega,MatSORType flag,double shift, 1021*bbf0e169SBarry Smith int its,Vec x) 1022*bbf0e169SBarry Smith { 1023*bbf0e169SBarry Smith int ierr; 1024*bbf0e169SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 1025*bbf0e169SBarry Smith PetscValidHeaderSpecific(b,VEC_COOKIE); PetscValidHeaderSpecific(x,VEC_COOKIE); 1026*bbf0e169SBarry Smith if (!mat->ops.relax) SETERRQ(PETSC_ERR_SUP,"MatRelax"); 1027*bbf0e169SBarry Smith if (!mat->assembled) SETERRQ(1,"MatRelax:Not for unassembled matrix"); 1028*bbf0e169SBarry Smith if (mat->factor) SETERRQ(1,"MatRelax:Not for factored matrix"); 1029*bbf0e169SBarry Smith if (mat->N != x->N) SETERRQ(PETSC_ERR_SIZ,"MatRelax:Mat mat,Vec x: global dim"); 1030*bbf0e169SBarry Smith if (mat->M != b->N) SETERRQ(PETSC_ERR_SIZ,"MatRelax:Mat mat,Vec b: global dim"); 1031*bbf0e169SBarry Smith if (mat->m != b->n) SETERRQ(PETSC_ERR_SIZ,"MatRelax:Mat mat,Vec b: local dim"); 1032*bbf0e169SBarry Smith 1033*bbf0e169SBarry Smith PLogEventBegin(MAT_Relax,mat,b,x,0); 1034*bbf0e169SBarry Smith ierr =(*mat->ops.relax)(mat,b,omega,flag,shift,its,x); CHKERRQ(ierr); 1035*bbf0e169SBarry Smith PLogEventEnd(MAT_Relax,mat,b,x,0); 1036*bbf0e169SBarry Smith return 0; 1037*bbf0e169SBarry Smith } 1038*bbf0e169SBarry Smith 1039*bbf0e169SBarry Smith /* 1040*bbf0e169SBarry Smith Default matrix copy routine. 1041*bbf0e169SBarry Smith */ 1042*bbf0e169SBarry Smith int MatCopy_Basic(Mat A,Mat B) 1043*bbf0e169SBarry Smith { 1044*bbf0e169SBarry Smith int ierr,i,rstart,rend,nz,*cwork; 1045*bbf0e169SBarry Smith Scalar *vwork; 1046*bbf0e169SBarry Smith 1047*bbf0e169SBarry Smith ierr = MatZeroEntries(B); CHKERRQ(ierr); 1048*bbf0e169SBarry Smith ierr = MatGetOwnershipRange(A,&rstart,&rend); CHKERRQ(ierr); 1049*bbf0e169SBarry Smith for (i=rstart; i<rend; i++) { 1050*bbf0e169SBarry Smith ierr = MatGetRow(A,i,&nz,&cwork,&vwork); CHKERRQ(ierr); 1051*bbf0e169SBarry Smith ierr = MatSetValues(B,1,&i,nz,cwork,vwork,INSERT_VALUES); CHKERRQ(ierr); 1052*bbf0e169SBarry Smith ierr = MatRestoreRow(A,i,&nz,&cwork,&vwork); CHKERRQ(ierr); 1053*bbf0e169SBarry Smith } 1054*bbf0e169SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 1055*bbf0e169SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 1056*bbf0e169SBarry Smith return 0; 1057*bbf0e169SBarry Smith } 1058*bbf0e169SBarry Smith 1059*bbf0e169SBarry Smith /*@C 1060*bbf0e169SBarry Smith MatCopy - Copys a matrix to another matrix. 1061*bbf0e169SBarry Smith 1062*bbf0e169SBarry Smith Input Parameters: 1063*bbf0e169SBarry Smith . A - the matrix 1064*bbf0e169SBarry Smith 1065*bbf0e169SBarry Smith Output Parameter: 1066*bbf0e169SBarry Smith . B - where the copy is put 1067*bbf0e169SBarry Smith 1068*bbf0e169SBarry Smith Notes: 1069*bbf0e169SBarry Smith MatCopy() copies the matrix entries of a matrix to another existing 1070*bbf0e169SBarry Smith matrix (after first zeroing the second matrix). A related routine is 1071*bbf0e169SBarry Smith MatConvert(), which first creates a new matrix and then copies the data. 1072*bbf0e169SBarry Smith 1073*bbf0e169SBarry Smith .keywords: matrix, copy, convert 1074*bbf0e169SBarry Smith 1075*bbf0e169SBarry Smith .seealso: MatConvert() 1076*bbf0e169SBarry Smith @*/ 1077*bbf0e169SBarry Smith int MatCopy(Mat A,Mat B) 1078*bbf0e169SBarry Smith { 1079*bbf0e169SBarry Smith int ierr; 1080*bbf0e169SBarry Smith PetscValidHeaderSpecific(A,MAT_COOKIE); PetscValidHeaderSpecific(B,MAT_COOKIE); 1081*bbf0e169SBarry Smith if (!A->assembled) SETERRQ(1,"MatCopy:Not for unassembled matrix"); 1082*bbf0e169SBarry Smith if (A->factor) SETERRQ(1,"MatCopy:Not for factored matrix"); 1083*bbf0e169SBarry Smith if (A->M != B->M || A->N != B->N) SETERRQ(PETSC_ERR_SIZ,"MatCopy:Mat A,Mat B: global dim"); 1084*bbf0e169SBarry Smith 1085*bbf0e169SBarry Smith PLogEventBegin(MAT_Copy,A,B,0,0); 1086*bbf0e169SBarry Smith if (A->ops.copy) { 1087*bbf0e169SBarry Smith ierr = (*A->ops.copy)(A,B); CHKERRQ(ierr); 1088*bbf0e169SBarry Smith } 1089*bbf0e169SBarry Smith else { /* generic conversion */ 1090*bbf0e169SBarry Smith ierr = MatCopy_Basic(A,B); CHKERRQ(ierr); 1091*bbf0e169SBarry Smith } 1092*bbf0e169SBarry Smith PLogEventEnd(MAT_Copy,A,B,0,0); 1093*bbf0e169SBarry Smith return 0; 1094*bbf0e169SBarry Smith } 1095*bbf0e169SBarry Smith 1096*bbf0e169SBarry Smith /*@C 1097*bbf0e169SBarry Smith MatConvert - Converts a matrix to another matrix, either of the same 1098*bbf0e169SBarry Smith or different type. 1099*bbf0e169SBarry Smith 1100*bbf0e169SBarry Smith Input Parameters: 1101*bbf0e169SBarry Smith . mat - the matrix 1102*bbf0e169SBarry Smith . newtype - new matrix type. Use MATSAME to create a new matrix of the 1103*bbf0e169SBarry Smith same type as the original matrix. 1104*bbf0e169SBarry Smith 1105*bbf0e169SBarry Smith Output Parameter: 1106*bbf0e169SBarry Smith . M - pointer to place new matrix 1107*bbf0e169SBarry Smith 1108*bbf0e169SBarry Smith Notes: 1109*bbf0e169SBarry Smith MatConvert() first creates a new matrix and then copies the data from 1110*bbf0e169SBarry Smith the first matrix. A related routine is MatCopy(), which copies the matrix 1111*bbf0e169SBarry Smith entries of one matrix to another already existing matrix context. 1112*bbf0e169SBarry Smith 1113*bbf0e169SBarry Smith .keywords: matrix, copy, convert 1114*bbf0e169SBarry Smith 1115*bbf0e169SBarry Smith .seealso: MatCopy() 1116*bbf0e169SBarry Smith @*/ 1117*bbf0e169SBarry Smith int MatConvert(Mat mat,MatType newtype,Mat *M) 1118*bbf0e169SBarry Smith { 1119*bbf0e169SBarry Smith int ierr; 1120*bbf0e169SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 1121*bbf0e169SBarry Smith if (!M) SETERRQ(1,"MatConvert:Bad new matrix address"); 1122*bbf0e169SBarry Smith if (!mat->assembled) SETERRQ(1,"MatConvert:Not for unassembled matrix"); 1123*bbf0e169SBarry Smith if (mat->factor) SETERRQ(1,"MatConvert:Not for factored matrix"); 1124*bbf0e169SBarry Smith 1125*bbf0e169SBarry Smith PLogEventBegin(MAT_Convert,mat,0,0,0); 1126*bbf0e169SBarry Smith if (newtype == mat->type || newtype == MATSAME) { 1127*bbf0e169SBarry Smith if (mat->ops.convertsametype) { /* customized copy */ 1128*bbf0e169SBarry Smith ierr = (*mat->ops.convertsametype)(mat,M,COPY_VALUES); CHKERRQ(ierr); 1129*bbf0e169SBarry Smith } 1130*bbf0e169SBarry Smith else { /* generic conversion */ 1131*bbf0e169SBarry Smith ierr = MatConvert_Basic(mat,newtype,M); CHKERRQ(ierr); 1132*bbf0e169SBarry Smith } 1133*bbf0e169SBarry Smith } 1134*bbf0e169SBarry Smith else if (mat->ops.convert) { /* customized conversion */ 1135*bbf0e169SBarry Smith ierr = (*mat->ops.convert)(mat,newtype,M); CHKERRQ(ierr); 1136*bbf0e169SBarry Smith } 1137*bbf0e169SBarry Smith else { /* generic conversion */ 1138*bbf0e169SBarry Smith ierr = MatConvert_Basic(mat,newtype,M); CHKERRQ(ierr); 1139*bbf0e169SBarry Smith } 1140*bbf0e169SBarry Smith PLogEventEnd(MAT_Convert,mat,0,0,0); 1141*bbf0e169SBarry Smith return 0; 1142*bbf0e169SBarry Smith } 1143*bbf0e169SBarry Smith 1144*bbf0e169SBarry Smith /*@ 1145*bbf0e169SBarry Smith MatGetDiagonal - Gets the diagonal of a matrix. 1146*bbf0e169SBarry Smith 1147*bbf0e169SBarry Smith Input Parameters: 1148*bbf0e169SBarry Smith . mat - the matrix 1149*bbf0e169SBarry Smith . v - the vector for storing the diagonal 1150*bbf0e169SBarry Smith 1151*bbf0e169SBarry Smith Output Parameter: 1152*bbf0e169SBarry Smith . v - the diagonal of the matrix 1153*bbf0e169SBarry Smith 1154*bbf0e169SBarry Smith .keywords: matrix, get, diagonal 1155*bbf0e169SBarry Smith @*/ 1156*bbf0e169SBarry Smith int MatGetDiagonal(Mat mat,Vec v) 1157*bbf0e169SBarry Smith { 1158*bbf0e169SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE);PetscValidHeaderSpecific(v,VEC_COOKIE); 1159*bbf0e169SBarry Smith if (!mat->assembled) SETERRQ(1,"MatGetDiagonal:Not for unassembled matrix"); 1160*bbf0e169SBarry Smith if (mat->factor) SETERRQ(1,"MatGetDiagonal:Not for factored matrix"); 1161*bbf0e169SBarry Smith if (mat->ops.getdiagonal) return (*mat->ops.getdiagonal)(mat,v); 1162*bbf0e169SBarry Smith SETERRQ(PETSC_ERR_SUP,"MatGetDiagonal"); 1163*bbf0e169SBarry Smith } 1164*bbf0e169SBarry Smith 1165*bbf0e169SBarry Smith /*@C 1166*bbf0e169SBarry Smith MatTranspose - Computes an in-place or out-of-place transpose of a matrix. 1167*bbf0e169SBarry Smith 1168*bbf0e169SBarry Smith Input Parameter: 1169*bbf0e169SBarry Smith . mat - the matrix to transpose 1170*bbf0e169SBarry Smith 1171*bbf0e169SBarry Smith Output Parameters: 1172*bbf0e169SBarry Smith . B - the transpose (or pass in PETSC_NULL for an in-place transpose) 1173*bbf0e169SBarry Smith 1174*bbf0e169SBarry Smith .keywords: matrix, transpose 1175*bbf0e169SBarry Smith 1176*bbf0e169SBarry Smith .seealso: MatMultTrans(), MatMultTransAdd() 1177*bbf0e169SBarry Smith @*/ 1178*bbf0e169SBarry Smith int MatTranspose(Mat mat,Mat *B) 1179*bbf0e169SBarry Smith { 1180*bbf0e169SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 1181*bbf0e169SBarry Smith if (!mat->assembled) SETERRQ(1,"MatTranspose:Not for unassembled matrix"); 1182*bbf0e169SBarry Smith if (mat->factor) SETERRQ(1,"MatTranspose:Not for factored matrix"); 1183*bbf0e169SBarry Smith if (mat->ops.transpose) return (*mat->ops.transpose)(mat,B); 1184*bbf0e169SBarry Smith SETERRQ(PETSC_ERR_SUP,"MatTranspose"); 1185*bbf0e169SBarry Smith } 1186*bbf0e169SBarry Smith 1187*bbf0e169SBarry Smith /*@ 1188*bbf0e169SBarry Smith MatEqual - Compares two matrices. 1189*bbf0e169SBarry Smith 1190*bbf0e169SBarry Smith Input Parameters: 1191*bbf0e169SBarry Smith . A - the first matrix 1192*bbf0e169SBarry Smith . B - the second matrix 1193*bbf0e169SBarry Smith 1194*bbf0e169SBarry Smith Output Parameter: 1195*bbf0e169SBarry Smith . flg : PETSC_TRUE if the matrices are equal; 1196*bbf0e169SBarry Smith PETSC_FALSE otherwise. 1197*bbf0e169SBarry Smith 1198*bbf0e169SBarry Smith .keywords: matrix, equal, equivalent 1199*bbf0e169SBarry Smith @*/ 1200*bbf0e169SBarry Smith int MatEqual(Mat A,Mat B,PetscTruth *flg) 1201*bbf0e169SBarry Smith { 1202*bbf0e169SBarry Smith PetscValidHeaderSpecific(A,MAT_COOKIE); PetscValidHeaderSpecific(B,MAT_COOKIE); 1203*bbf0e169SBarry Smith PetscValidIntPointer(flg); 1204*bbf0e169SBarry Smith if (!A->assembled) SETERRQ(1,"MatEqual:Not for unassembled matrix"); 1205*bbf0e169SBarry Smith if (!B->assembled) SETERRQ(1,"MatEqual:Not for unassembled matrix"); 1206*bbf0e169SBarry Smith if (A->M != B->M || A->N != B->N) SETERRQ(PETSC_ERR_SIZ,"MatCopy:Mat A,Mat B: global dim"); 1207*bbf0e169SBarry Smith if (A->ops.equal) return (*A->ops.equal)(A,B,flg); 1208*bbf0e169SBarry Smith SETERRQ(PETSC_ERR_SUP,"MatEqual"); 1209*bbf0e169SBarry Smith } 1210*bbf0e169SBarry Smith 1211*bbf0e169SBarry Smith /*@ 1212*bbf0e169SBarry Smith MatDiagonalScale - Scales a matrix on the left and right by diagonal 1213*bbf0e169SBarry Smith matrices that are stored as vectors. Either of the two scaling 1214*bbf0e169SBarry Smith matrices can be PETSC_NULL. 1215*bbf0e169SBarry Smith 1216*bbf0e169SBarry Smith Input Parameters: 1217*bbf0e169SBarry Smith . mat - the matrix to be scaled 1218*bbf0e169SBarry Smith . l - the left scaling vector (or PETSC_NULL) 1219*bbf0e169SBarry Smith . r - the right scaling vector (or PETSC_NULL) 1220*bbf0e169SBarry Smith 1221*bbf0e169SBarry Smith Notes: 1222*bbf0e169SBarry Smith MatDiagonalScale() computes A <- LAR, where 1223*bbf0e169SBarry Smith $ L = a diagonal matrix 1224*bbf0e169SBarry Smith $ R = a diagonal matrix 1225*bbf0e169SBarry Smith 1226*bbf0e169SBarry Smith .keywords: matrix, diagonal, scale 1227*bbf0e169SBarry Smith 1228*bbf0e169SBarry Smith .seealso: MatDiagonalScale() 1229*bbf0e169SBarry Smith @*/ 1230*bbf0e169SBarry Smith int MatDiagonalScale(Mat mat,Vec l,Vec r) 1231*bbf0e169SBarry Smith { 1232*bbf0e169SBarry Smith int ierr; 1233*bbf0e169SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 1234*bbf0e169SBarry Smith if (!mat->ops.diagonalscale) SETERRQ(PETSC_ERR_SUP,"MatDiagonalScale"); 1235*bbf0e169SBarry Smith if (l) PetscValidHeaderSpecific(l,VEC_COOKIE); 1236*bbf0e169SBarry Smith if (r) PetscValidHeaderSpecific(r,VEC_COOKIE); 1237*bbf0e169SBarry Smith if (!mat->assembled) SETERRQ(1,"MatDiagonalScale:Not for unassembled matrix"); 1238*bbf0e169SBarry Smith if (mat->factor) SETERRQ(1,"MatDiagonalScale:Not for factored matrix"); 1239*bbf0e169SBarry Smith 1240*bbf0e169SBarry Smith PLogEventBegin(MAT_Scale,mat,0,0,0); 1241*bbf0e169SBarry Smith ierr = (*mat->ops.diagonalscale)(mat,l,r); CHKERRQ(ierr); 1242*bbf0e169SBarry Smith PLogEventEnd(MAT_Scale,mat,0,0,0); 1243*bbf0e169SBarry Smith return 0; 1244*bbf0e169SBarry Smith } 1245*bbf0e169SBarry Smith 1246*bbf0e169SBarry Smith /*@ 1247*bbf0e169SBarry Smith MatScale - Scales all elements of a matrix by a given number. 1248*bbf0e169SBarry Smith 1249*bbf0e169SBarry Smith Input Parameters: 1250*bbf0e169SBarry Smith . mat - the matrix to be scaled 1251*bbf0e169SBarry Smith . a - the scaling value 1252*bbf0e169SBarry Smith 1253*bbf0e169SBarry Smith Output Parameter: 1254*bbf0e169SBarry Smith . mat - the scaled matrix 1255*bbf0e169SBarry Smith 1256*bbf0e169SBarry Smith .keywords: matrix, scale 1257*bbf0e169SBarry Smith 1258*bbf0e169SBarry Smith .seealso: MatDiagonalScale() 1259*bbf0e169SBarry Smith @*/ 1260*bbf0e169SBarry Smith int MatScale(Scalar *a,Mat mat) 1261*bbf0e169SBarry Smith { 1262*bbf0e169SBarry Smith int ierr; 1263*bbf0e169SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 1264*bbf0e169SBarry Smith PetscValidScalarPointer(a); 1265*bbf0e169SBarry Smith if (!mat->ops.scale) SETERRQ(PETSC_ERR_SUP,"MatScale"); 1266*bbf0e169SBarry Smith if (!mat->assembled) SETERRQ(1,"MatScale:Not for unassembled matrix"); 1267*bbf0e169SBarry Smith if (mat->factor) SETERRQ(1,"MatScale:Not for factored matrix"); 1268*bbf0e169SBarry Smith 1269*bbf0e169SBarry Smith PLogEventBegin(MAT_Scale,mat,0,0,0); 1270*bbf0e169SBarry Smith ierr = (*mat->ops.scale)(a,mat); CHKERRQ(ierr); 1271*bbf0e169SBarry Smith PLogEventEnd(MAT_Scale,mat,0,0,0); 1272*bbf0e169SBarry Smith return 0; 1273*bbf0e169SBarry Smith } 1274*bbf0e169SBarry Smith 1275*bbf0e169SBarry Smith /*@ 1276*bbf0e169SBarry Smith MatNorm - Calculates various norms of a matrix. 1277*bbf0e169SBarry Smith 1278*bbf0e169SBarry Smith Input Parameters: 1279*bbf0e169SBarry Smith . mat - the matrix 1280*bbf0e169SBarry Smith . type - the type of norm, NORM_1, NORM_2, NORM_FROBENIUS, NORM_INFINITY 1281*bbf0e169SBarry Smith 1282*bbf0e169SBarry Smith Output Parameters: 1283*bbf0e169SBarry Smith . norm - the resulting norm 1284*bbf0e169SBarry Smith 1285*bbf0e169SBarry Smith .keywords: matrix, norm, Frobenius 1286*bbf0e169SBarry Smith @*/ 1287*bbf0e169SBarry Smith int MatNorm(Mat mat,NormType type,double *norm) 1288*bbf0e169SBarry Smith { 1289*bbf0e169SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 1290*bbf0e169SBarry Smith PetscValidScalarPointer(norm); 1291*bbf0e169SBarry Smith 1292*bbf0e169SBarry Smith if (!norm) SETERRQ(1,"MatNorm:bad addess for value"); 1293*bbf0e169SBarry Smith if (!mat->assembled) SETERRQ(1,"MatNorm:Not for unassembled matrix"); 1294*bbf0e169SBarry Smith if (mat->factor) SETERRQ(1,"MatNorm:Not for factored matrix"); 1295*bbf0e169SBarry Smith if (mat->ops.norm) return (*mat->ops.norm)(mat,type,norm); 1296*bbf0e169SBarry Smith SETERRQ(PETSC_ERR_SUP,"MatNorm:Not for this matrix type"); 1297*bbf0e169SBarry Smith } 1298*bbf0e169SBarry Smith 1299*bbf0e169SBarry Smith /*@ 1300*bbf0e169SBarry Smith MatAssemblyBegin - Begins assembling the matrix. This routine should 1301*bbf0e169SBarry Smith be called after completing all calls to MatSetValues(). 1302*bbf0e169SBarry Smith 1303*bbf0e169SBarry Smith Input Parameters: 1304*bbf0e169SBarry Smith . mat - the matrix 1305*bbf0e169SBarry Smith . type - type of assembly, either MAT_FLUSH_ASSEMBLY or MAT_FINAL_ASSEMBLY 1306*bbf0e169SBarry Smith 1307*bbf0e169SBarry Smith Notes: 1308*bbf0e169SBarry Smith MatSetValues() generally caches the values. The matrix is ready to 1309*bbf0e169SBarry Smith use only after MatAssemblyBegin() and MatAssemblyEnd() have been called. 1310*bbf0e169SBarry Smith Use MAT_FLUSH_ASSEMBLY when switching between ADD_VALUES and INSERT_VALUES 1311*bbf0e169SBarry Smith in MatSetValues(); use MAT_FINAL_ASSEMBLY for the final assembly before 1312*bbf0e169SBarry Smith using the matrix. 1313*bbf0e169SBarry Smith 1314*bbf0e169SBarry Smith .keywords: matrix, assembly, assemble, begin 1315*bbf0e169SBarry Smith 1316*bbf0e169SBarry Smith .seealso: MatAssemblyEnd(), MatSetValues() 1317*bbf0e169SBarry Smith @*/ 1318*bbf0e169SBarry Smith int MatAssemblyBegin(Mat mat,MatAssemblyType type) 1319*bbf0e169SBarry Smith { 1320*bbf0e169SBarry Smith int ierr; 1321*bbf0e169SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 1322*bbf0e169SBarry Smith if (mat->factor) SETERRQ(1,"MatAssemblyBegin:Not for factored matrix"); 1323*bbf0e169SBarry Smith if (mat->assembled) { 1324*bbf0e169SBarry Smith mat->was_assembled = PETSC_TRUE; 1325*bbf0e169SBarry Smith mat->assembled = PETSC_FALSE; 1326*bbf0e169SBarry Smith } 1327*bbf0e169SBarry Smith PLogEventBegin(MAT_AssemblyBegin,mat,0,0,0); 1328*bbf0e169SBarry Smith if (mat->ops.assemblybegin){ierr = (*mat->ops.assemblybegin)(mat,type);CHKERRQ(ierr);} 1329*bbf0e169SBarry Smith PLogEventEnd(MAT_AssemblyBegin,mat,0,0,0); 1330*bbf0e169SBarry Smith return 0; 1331*bbf0e169SBarry Smith } 1332*bbf0e169SBarry Smith 1333*bbf0e169SBarry Smith /*@ 1334*bbf0e169SBarry Smith MatAssemblyEnd - Completes assembling the matrix. This routine should 1335*bbf0e169SBarry Smith be called after MatAssemblyBegin(). 1336*bbf0e169SBarry Smith 1337*bbf0e169SBarry Smith Input Parameters: 1338*bbf0e169SBarry Smith . mat - the matrix 1339*bbf0e169SBarry Smith . type - type of assembly, either MAT_FLUSH_ASSEMBLY or MAT_FINAL_ASSEMBLY 1340*bbf0e169SBarry Smith 1341*bbf0e169SBarry Smith Options Database Keys: 1342*bbf0e169SBarry Smith $ -mat_view_info : Prints info on matrix at 1343*bbf0e169SBarry Smith $ conclusion of MatEndAssembly() 1344*bbf0e169SBarry Smith $ -mat_view_info_detailed: Prints more detailed info. 1345*bbf0e169SBarry Smith $ -mat_view : Prints matrix in ASCII format. 1346*bbf0e169SBarry Smith $ -mat_view_matlab : Prints matrix in Matlab format. 1347*bbf0e169SBarry Smith $ -mat_view_draw : Draws nonzero structure of matrix, 1348*bbf0e169SBarry Smith $ using MatView() and DrawOpenX(). 1349*bbf0e169SBarry Smith $ -display <name> : Set display name (default is host) 1350*bbf0e169SBarry Smith $ -draw_pause <sec> : Set number of seconds to pause after display 1351*bbf0e169SBarry Smith 1352*bbf0e169SBarry Smith Notes: 1353*bbf0e169SBarry Smith MatSetValues() generally caches the values. The matrix is ready to 1354*bbf0e169SBarry Smith use only after MatAssemblyBegin() and MatAssemblyEnd() have been called. 1355*bbf0e169SBarry Smith Use MAT_FLUSH_ASSEMBLY when switching between ADD_VALUES and INSERT_VALUES 1356*bbf0e169SBarry Smith in MatSetValues(); use MAT_FINAL_ASSEMBLY for the final assembly before 1357*bbf0e169SBarry Smith using the matrix. 1358*bbf0e169SBarry Smith 1359*bbf0e169SBarry Smith .keywords: matrix, assembly, assemble, end 1360*bbf0e169SBarry Smith 1361*bbf0e169SBarry Smith .seealso: MatAssemblyBegin(), MatSetValues(), DrawOpenX(), MatView() 1362*bbf0e169SBarry Smith @*/ 1363*bbf0e169SBarry Smith int MatAssemblyEnd(Mat mat,MatAssemblyType type) 1364*bbf0e169SBarry Smith { 1365*bbf0e169SBarry Smith int ierr,flg; 1366*bbf0e169SBarry Smith static int inassm = 0; 1367*bbf0e169SBarry Smith 1368*bbf0e169SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 1369*bbf0e169SBarry Smith inassm++; 1370*bbf0e169SBarry Smith PLogEventBegin(MAT_AssemblyEnd,mat,0,0,0); 1371*bbf0e169SBarry Smith if (mat->ops.assemblyend) { 1372*bbf0e169SBarry Smith ierr = (*mat->ops.assemblyend)(mat,type); CHKERRQ(ierr); 1373*bbf0e169SBarry Smith } 1374*bbf0e169SBarry Smith mat->assembled = PETSC_TRUE; mat->num_ass++; 1375*bbf0e169SBarry Smith PLogEventEnd(MAT_AssemblyEnd,mat,0,0,0); 1376*bbf0e169SBarry Smith 1377*bbf0e169SBarry Smith if (inassm == 1) { 1378*bbf0e169SBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_view_info",&flg); CHKERRQ(ierr); 1379*bbf0e169SBarry Smith if (flg) { 1380*bbf0e169SBarry Smith Viewer viewer; 1381*bbf0e169SBarry Smith ierr = ViewerFileOpenASCII(mat->comm,"stdout",&viewer);CHKERRQ(ierr); 1382*bbf0e169SBarry Smith ierr = ViewerSetFormat(viewer,VIEWER_FORMAT_ASCII_INFO,0);CHKERRQ(ierr); 1383*bbf0e169SBarry Smith ierr = MatView(mat,viewer); CHKERRQ(ierr); 1384*bbf0e169SBarry Smith ierr = ViewerDestroy(viewer); CHKERRQ(ierr); 1385*bbf0e169SBarry Smith } 1386*bbf0e169SBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_view_info_detailed",&flg);CHKERRQ(ierr); 1387*bbf0e169SBarry Smith if (flg) { 1388*bbf0e169SBarry Smith Viewer viewer; 1389*bbf0e169SBarry Smith ierr = ViewerFileOpenASCII(mat->comm,"stdout",&viewer);CHKERRQ(ierr); 1390*bbf0e169SBarry Smith ierr = ViewerSetFormat(viewer,VIEWER_FORMAT_ASCII_INFO_LONG,0);CHKERRQ(ierr); 1391*bbf0e169SBarry Smith ierr = MatView(mat,viewer); CHKERRQ(ierr); 1392*bbf0e169SBarry Smith ierr = ViewerDestroy(viewer); CHKERRQ(ierr); 1393*bbf0e169SBarry Smith } 1394*bbf0e169SBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_view",&flg); CHKERRQ(ierr); 1395*bbf0e169SBarry Smith if (flg) { 1396*bbf0e169SBarry Smith Viewer viewer; 1397*bbf0e169SBarry Smith ierr = ViewerFileOpenASCII(mat->comm,"stdout",&viewer);CHKERRQ(ierr); 1398*bbf0e169SBarry Smith ierr = MatView(mat,viewer); CHKERRQ(ierr); 1399*bbf0e169SBarry Smith ierr = ViewerDestroy(viewer); CHKERRQ(ierr); 1400*bbf0e169SBarry Smith } 1401*bbf0e169SBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_view_matlab",&flg); CHKERRQ(ierr); 1402*bbf0e169SBarry Smith if (flg) { 1403*bbf0e169SBarry Smith Viewer viewer; 1404*bbf0e169SBarry Smith ierr = ViewerFileOpenASCII(mat->comm,"stdout",&viewer);CHKERRQ(ierr); 1405*bbf0e169SBarry Smith ierr = ViewerSetFormat(viewer,VIEWER_FORMAT_ASCII_MATLAB,"M");CHKERRQ(ierr); 1406*bbf0e169SBarry Smith ierr = MatView(mat,viewer); CHKERRQ(ierr); 1407*bbf0e169SBarry Smith ierr = ViewerDestroy(viewer); CHKERRQ(ierr); 1408*bbf0e169SBarry Smith } 1409*bbf0e169SBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_view_draw",&flg); CHKERRQ(ierr); 1410*bbf0e169SBarry Smith if (flg) { 1411*bbf0e169SBarry Smith Viewer viewer; 1412*bbf0e169SBarry Smith ierr = ViewerDrawOpenX(mat->comm,0,0,0,0,300,300,&viewer); CHKERRQ(ierr); 1413*bbf0e169SBarry Smith ierr = MatView(mat,viewer); CHKERRQ(ierr); 1414*bbf0e169SBarry Smith ierr = ViewerFlush(viewer); CHKERRQ(ierr); 1415*bbf0e169SBarry Smith ierr = ViewerDestroy(viewer); CHKERRQ(ierr); 1416*bbf0e169SBarry Smith } 1417*bbf0e169SBarry Smith } 1418*bbf0e169SBarry Smith inassm--; 1419*bbf0e169SBarry Smith return 0; 1420*bbf0e169SBarry Smith } 1421*bbf0e169SBarry Smith 1422*bbf0e169SBarry Smith /*@ 1423*bbf0e169SBarry Smith MatCompress - Tries to store the matrix in as little space as 1424*bbf0e169SBarry Smith possible. May fail if memory is already fully used, since it 1425*bbf0e169SBarry Smith tries to allocate new space. 1426*bbf0e169SBarry Smith 1427*bbf0e169SBarry Smith Input Parameters: 1428*bbf0e169SBarry Smith . mat - the matrix 1429*bbf0e169SBarry Smith 1430*bbf0e169SBarry Smith .keywords: matrix, compress 1431*bbf0e169SBarry Smith @*/ 1432*bbf0e169SBarry Smith int MatCompress(Mat mat) 1433*bbf0e169SBarry Smith { 1434*bbf0e169SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 1435*bbf0e169SBarry Smith if (mat->ops.compress) return (*mat->ops.compress)(mat); 1436*bbf0e169SBarry Smith return 0; 1437*bbf0e169SBarry Smith } 1438*bbf0e169SBarry Smith /*@ 1439*bbf0e169SBarry Smith MatSetOption - Sets a parameter option for a matrix. Some options 1440*bbf0e169SBarry Smith may be specific to certain storage formats. Some options 1441*bbf0e169SBarry Smith determine how values will be inserted (or added). Sorted, 1442*bbf0e169SBarry Smith row-oriented input will generally assemble the fastest. The default 1443*bbf0e169SBarry Smith is row-oriented, nonsorted input. 1444*bbf0e169SBarry Smith 1445*bbf0e169SBarry Smith Input Parameters: 1446*bbf0e169SBarry Smith . mat - the matrix 1447*bbf0e169SBarry Smith . option - the option, one of the following: 1448*bbf0e169SBarry Smith $ MAT_ROW_ORIENTED 1449*bbf0e169SBarry Smith $ MAT_COLUMN_ORIENTED, 1450*bbf0e169SBarry Smith $ MAT_ROWS_SORTED, 1451*bbf0e169SBarry Smith $ MAT_COLUMNS_SORTED, 1452*bbf0e169SBarry Smith $ MAT_NO_NEW_NONZERO_LOCATIONS, 1453*bbf0e169SBarry Smith $ MAT_YES_NEW_NONZERO_LOCATIONS, 1454*bbf0e169SBarry Smith $ MAT_SYMMETRIC, 1455*bbf0e169SBarry Smith $ MAT_STRUCTURALLY_SYMMETRIC, 1456*bbf0e169SBarry Smith $ MAT_NO_NEW_DIAGONALS, 1457*bbf0e169SBarry Smith $ MAT_YES_NEW_DIAGONALS, 1458*bbf0e169SBarry Smith $ and possibly others. 1459*bbf0e169SBarry Smith 1460*bbf0e169SBarry Smith Notes: 1461*bbf0e169SBarry Smith Some options are relevant only for particular matrix types and 1462*bbf0e169SBarry Smith are thus ignored by others. Other options are not supported by 1463*bbf0e169SBarry Smith certain matrix types and will generate an error message if set. 1464*bbf0e169SBarry Smith 1465*bbf0e169SBarry Smith If using a Fortran 77 module to compute a matrix, one may need to 1466*bbf0e169SBarry Smith use the column-oriented option (or convert to the row-oriented 1467*bbf0e169SBarry Smith format). 1468*bbf0e169SBarry Smith 1469*bbf0e169SBarry Smith MAT_NO_NEW_NONZERO_LOCATIONS indicates that any add or insertion 1470*bbf0e169SBarry Smith that will generate a new entry in the nonzero structure is ignored. 1471*bbf0e169SBarry Smith What this means is if memory is not allocated for this particular 1472*bbf0e169SBarry Smith lot, then the insertion is ignored. For dense matrices, where 1473*bbf0e169SBarry Smith the entire array is allocated, no entries are ever ignored. 1474*bbf0e169SBarry Smith 1475*bbf0e169SBarry Smith .keywords: matrix, option, row-oriented, column-oriented, sorted, nonzero 1476*bbf0e169SBarry Smith @*/ 1477*bbf0e169SBarry Smith int MatSetOption(Mat mat,MatOption op) 1478*bbf0e169SBarry Smith { 1479*bbf0e169SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 1480*bbf0e169SBarry Smith if (mat->ops.setoption) return (*mat->ops.setoption)(mat,op); 1481*bbf0e169SBarry Smith return 0; 1482*bbf0e169SBarry Smith } 1483*bbf0e169SBarry Smith 1484*bbf0e169SBarry Smith /*@ 1485*bbf0e169SBarry Smith MatZeroEntries - Zeros all entries of a matrix. For sparse matrices 1486*bbf0e169SBarry Smith this routine retains the old nonzero structure. 1487*bbf0e169SBarry Smith 1488*bbf0e169SBarry Smith Input Parameters: 1489*bbf0e169SBarry Smith . mat - the matrix 1490*bbf0e169SBarry Smith 1491*bbf0e169SBarry Smith .keywords: matrix, zero, entries 1492*bbf0e169SBarry Smith 1493*bbf0e169SBarry Smith .seealso: MatZeroRows() 1494*bbf0e169SBarry Smith @*/ 1495*bbf0e169SBarry Smith int MatZeroEntries(Mat mat) 1496*bbf0e169SBarry Smith { 1497*bbf0e169SBarry Smith int ierr; 1498*bbf0e169SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 1499*bbf0e169SBarry Smith if (mat->factor) SETERRQ(1,"MatZeroEntries:Not for factored matrix"); 1500*bbf0e169SBarry Smith if (!mat->ops.zeroentries) SETERRQ(PETSC_ERR_SUP,"MatZeroEntries"); 1501*bbf0e169SBarry Smith 1502*bbf0e169SBarry Smith PLogEventBegin(MAT_ZeroEntries,mat,0,0,0); 1503*bbf0e169SBarry Smith ierr = (*mat->ops.zeroentries)(mat); CHKERRQ(ierr); 1504*bbf0e169SBarry Smith PLogEventEnd(MAT_ZeroEntries,mat,0,0,0); 1505*bbf0e169SBarry Smith return 0; 1506*bbf0e169SBarry Smith } 1507*bbf0e169SBarry Smith 1508*bbf0e169SBarry Smith /*@ 1509*bbf0e169SBarry Smith MatZeroRows - Zeros all entries (except possibly the main diagonal) 1510*bbf0e169SBarry Smith of a set of rows of a matrix. 1511*bbf0e169SBarry Smith 1512*bbf0e169SBarry Smith Input Parameters: 1513*bbf0e169SBarry Smith . mat - the matrix 1514*bbf0e169SBarry Smith . is - index set of rows to remove 1515*bbf0e169SBarry Smith . diag - pointer to value put in all diagonals of eliminated rows. 1516*bbf0e169SBarry Smith Note that diag is not a pointer to an array, but merely a 1517*bbf0e169SBarry Smith pointer to a single value. 1518*bbf0e169SBarry Smith 1519*bbf0e169SBarry Smith Notes: 1520*bbf0e169SBarry Smith For the AIJ matrix formats this removes the old nonzero structure, 1521*bbf0e169SBarry Smith but does not release memory. For the dense and block diagonal 1522*bbf0e169SBarry Smith formats this does not alter the nonzero structure. 1523*bbf0e169SBarry Smith 1524*bbf0e169SBarry Smith The user can set a value in the diagonal entry (or for the AIJ and 1525*bbf0e169SBarry Smith row formats can optionally remove the main diagonal entry from the 1526*bbf0e169SBarry Smith nonzero structure as well, by passing a null pointer as the final 1527*bbf0e169SBarry Smith argument). 1528*bbf0e169SBarry Smith 1529*bbf0e169SBarry Smith .keywords: matrix, zero, rows, boundary conditions 1530*bbf0e169SBarry Smith 1531*bbf0e169SBarry Smith .seealso: MatZeroEntries(), 1532*bbf0e169SBarry Smith @*/ 1533*bbf0e169SBarry Smith int MatZeroRows(Mat mat,IS is, Scalar *diag) 1534*bbf0e169SBarry Smith { 1535*bbf0e169SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 1536*bbf0e169SBarry Smith if (!mat->assembled) SETERRQ(1,"MatZeroRows:Not for unassembled matrix"); 1537*bbf0e169SBarry Smith if (mat->factor) SETERRQ(1,"MatZeroRows:Not for factored matrix"); 1538*bbf0e169SBarry Smith if (mat->ops.zerorows) return (*mat->ops.zerorows)(mat,is,diag); 1539*bbf0e169SBarry Smith SETERRQ(PETSC_ERR_SUP,"MatZeroRows"); 1540*bbf0e169SBarry Smith } 1541*bbf0e169SBarry Smith 1542*bbf0e169SBarry Smith /*@ 1543*bbf0e169SBarry Smith MatGetSize - Returns the numbers of rows and columns in a matrix. 1544*bbf0e169SBarry Smith 1545*bbf0e169SBarry Smith Input Parameter: 1546*bbf0e169SBarry Smith . mat - the matrix 1547*bbf0e169SBarry Smith 1548*bbf0e169SBarry Smith Output Parameters: 1549*bbf0e169SBarry Smith . m - the number of global rows 1550*bbf0e169SBarry Smith . n - the number of global columns 1551*bbf0e169SBarry Smith 1552*bbf0e169SBarry Smith .keywords: matrix, dimension, size, rows, columns, global, get 1553*bbf0e169SBarry Smith 1554*bbf0e169SBarry Smith .seealso: MatGetLocalSize() 1555*bbf0e169SBarry Smith @*/ 1556*bbf0e169SBarry Smith int MatGetSize(Mat mat,int *m,int* n) 1557*bbf0e169SBarry Smith { 1558*bbf0e169SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 1559*bbf0e169SBarry Smith PetscValidIntPointer(m); 1560*bbf0e169SBarry Smith PetscValidIntPointer(n); 1561*bbf0e169SBarry Smith return (*mat->ops.getsize)(mat,m,n); 1562*bbf0e169SBarry Smith } 1563*bbf0e169SBarry Smith 1564*bbf0e169SBarry Smith /*@ 1565*bbf0e169SBarry Smith MatGetLocalSize - Returns the number of rows and columns in a matrix 1566*bbf0e169SBarry Smith stored locally. This information may be implementation dependent, so 1567*bbf0e169SBarry Smith use with care. 1568*bbf0e169SBarry Smith 1569*bbf0e169SBarry Smith Input Parameters: 1570*bbf0e169SBarry Smith . mat - the matrix 1571*bbf0e169SBarry Smith 1572*bbf0e169SBarry Smith Output Parameters: 1573*bbf0e169SBarry Smith . m - the number of local rows 1574*bbf0e169SBarry Smith . n - the number of local columns 1575*bbf0e169SBarry Smith 1576*bbf0e169SBarry Smith .keywords: matrix, dimension, size, local, rows, columns, get 1577*bbf0e169SBarry Smith 1578*bbf0e169SBarry Smith .seealso: MatGetSize() 1579*bbf0e169SBarry Smith @*/ 1580*bbf0e169SBarry Smith int MatGetLocalSize(Mat mat,int *m,int* n) 1581*bbf0e169SBarry Smith { 1582*bbf0e169SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 1583*bbf0e169SBarry Smith PetscValidIntPointer(m); 1584*bbf0e169SBarry Smith PetscValidIntPointer(n); 1585*bbf0e169SBarry Smith return (*mat->ops.getlocalsize)(mat,m,n); 1586*bbf0e169SBarry Smith } 1587*bbf0e169SBarry Smith 1588*bbf0e169SBarry Smith /*@ 1589*bbf0e169SBarry Smith MatGetOwnershipRange - Returns the range of matrix rows owned by 1590*bbf0e169SBarry Smith this processor, assuming that the matrix is laid out with the first 1591*bbf0e169SBarry Smith n1 rows on the first processor, the next n2 rows on the second, etc. 1592*bbf0e169SBarry Smith For certain parallel layouts this range may not be well defined. 1593*bbf0e169SBarry Smith 1594*bbf0e169SBarry Smith Input Parameters: 1595*bbf0e169SBarry Smith . mat - the matrix 1596*bbf0e169SBarry Smith 1597*bbf0e169SBarry Smith Output Parameters: 1598*bbf0e169SBarry Smith . m - the first local row 1599*bbf0e169SBarry Smith . n - one more then the last local row 1600*bbf0e169SBarry Smith 1601*bbf0e169SBarry Smith .keywords: matrix, get, range, ownership 1602*bbf0e169SBarry Smith @*/ 1603*bbf0e169SBarry Smith int MatGetOwnershipRange(Mat mat,int *m,int* n) 1604*bbf0e169SBarry Smith { 1605*bbf0e169SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 1606*bbf0e169SBarry Smith PetscValidIntPointer(m); 1607*bbf0e169SBarry Smith PetscValidIntPointer(n); 1608*bbf0e169SBarry Smith if (mat->ops.getownershiprange) return (*mat->ops.getownershiprange)(mat,m,n); 1609*bbf0e169SBarry Smith SETERRQ(PETSC_ERR_SUP,"MatGetOwnershipRange"); 1610*bbf0e169SBarry Smith } 1611*bbf0e169SBarry Smith 1612*bbf0e169SBarry Smith /*@ 1613*bbf0e169SBarry Smith MatILUFactorSymbolic - Performs symbolic ILU factorization of a matrix. 1614*bbf0e169SBarry Smith Uses levels of fill only, not drop tolerance. Use MatLUFactorNumeric() 1615*bbf0e169SBarry Smith to complete the factorization. 1616*bbf0e169SBarry Smith 1617*bbf0e169SBarry Smith Input Parameters: 1618*bbf0e169SBarry Smith . mat - the matrix 1619*bbf0e169SBarry Smith . row - row permutation 1620*bbf0e169SBarry Smith . column - column permutation 1621*bbf0e169SBarry Smith . fill - number of levels of fill 1622*bbf0e169SBarry Smith . f - expected fill as ratio of the original number of nonzeros, 1623*bbf0e169SBarry Smith for example 3.0; choosing this parameter well can result in 1624*bbf0e169SBarry Smith more efficient use of time and space. 1625*bbf0e169SBarry Smith 1626*bbf0e169SBarry Smith Output Parameters: 1627*bbf0e169SBarry Smith . fact - new matrix that has been symbolically factored 1628*bbf0e169SBarry Smith 1629*bbf0e169SBarry Smith Options Database Key: 1630*bbf0e169SBarry Smith $ -mat_ilu_fill <f>, where f is the fill ratio 1631*bbf0e169SBarry Smith 1632*bbf0e169SBarry Smith Notes: 1633*bbf0e169SBarry Smith See the file $(PETSC_DIR)/Performace for additional information about 1634*bbf0e169SBarry Smith choosing the fill factor for better efficiency. 1635*bbf0e169SBarry Smith 1636*bbf0e169SBarry Smith .keywords: matrix, factor, incomplete, ILU, symbolic, fill 1637*bbf0e169SBarry Smith 1638*bbf0e169SBarry Smith .seealso: MatLUFactorSymbolic(), MatLUFactorNumeric() 1639*bbf0e169SBarry Smith @*/ 1640*bbf0e169SBarry Smith int MatILUFactorSymbolic(Mat mat,IS row,IS col,double f,int fill,Mat *fact) 1641*bbf0e169SBarry Smith { 1642*bbf0e169SBarry Smith int ierr,flg; 1643*bbf0e169SBarry Smith 1644*bbf0e169SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 1645*bbf0e169SBarry Smith if (fill < 0) SETERRQ(1,"MatILUFactorSymbolic:Levels of fill negative"); 1646*bbf0e169SBarry Smith if (!fact) SETERRQ(1,"MatILUFactorSymbolic:Fact argument is missing"); 1647*bbf0e169SBarry Smith if (!mat->ops.ilufactorsymbolic) SETERRQ(PETSC_ERR_SUP,"MatILUFactorSymbolic"); 1648*bbf0e169SBarry Smith if (!mat->assembled) SETERRQ(1,"MatILUFactorSymbolic:Not for unassembled matrix"); 1649*bbf0e169SBarry Smith if (mat->factor) SETERRQ(1,"MatILUFactorSymbolic:Not for factored matrix"); 1650*bbf0e169SBarry Smith 1651*bbf0e169SBarry Smith ierr = OptionsGetDouble(PETSC_NULL,"-mat_ilu_fill",&f,&flg); CHKERRQ(ierr); 1652*bbf0e169SBarry Smith PLogEventBegin(MAT_ILUFactorSymbolic,mat,row,col,0); 1653*bbf0e169SBarry Smith ierr = (*mat->ops.ilufactorsymbolic)(mat,row,col,f,fill,fact); CHKERRQ(ierr); 1654*bbf0e169SBarry Smith PLogEventEnd(MAT_ILUFactorSymbolic,mat,row,col,0); 1655*bbf0e169SBarry Smith return 0; 1656*bbf0e169SBarry Smith } 1657*bbf0e169SBarry Smith 1658*bbf0e169SBarry Smith /*@ 1659*bbf0e169SBarry Smith MatIncompleteCholeskyFactorSymbolic - Performs symbolic incomplete 1660*bbf0e169SBarry Smith Cholesky factorization for a symmetric matrix. Use 1661*bbf0e169SBarry Smith MatCholeskyFactorNumeric() to complete the factorization. 1662*bbf0e169SBarry Smith 1663*bbf0e169SBarry Smith Input Parameters: 1664*bbf0e169SBarry Smith . mat - the matrix 1665*bbf0e169SBarry Smith . perm - row and column permutation 1666*bbf0e169SBarry Smith . fill - levels of fill 1667*bbf0e169SBarry Smith . f - expected fill as ratio of original fill 1668*bbf0e169SBarry Smith 1669*bbf0e169SBarry Smith Output Parameter: 1670*bbf0e169SBarry Smith . fact - the factored matrix 1671*bbf0e169SBarry Smith 1672*bbf0e169SBarry Smith Note: Currently only no-fill factorization is supported. 1673*bbf0e169SBarry Smith 1674*bbf0e169SBarry Smith .keywords: matrix, factor, incomplete, ICC, Cholesky, symbolic, fill 1675*bbf0e169SBarry Smith 1676*bbf0e169SBarry Smith .seealso: MatCholeskyFactorNumeric(), MatCholeskyFactor() 1677*bbf0e169SBarry Smith @*/ 1678*bbf0e169SBarry Smith int MatIncompleteCholeskyFactorSymbolic(Mat mat,IS perm,double f,int fill, 1679*bbf0e169SBarry Smith Mat *fact) 1680*bbf0e169SBarry Smith { 1681*bbf0e169SBarry Smith int ierr; 1682*bbf0e169SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 1683*bbf0e169SBarry Smith if (mat->factor) SETERRQ(1,"MatIncompleteCholeskyFactorSymbolic:Not for factored matrix"); 1684*bbf0e169SBarry Smith if (fill < 0) SETERRQ(1,"MatIncompleteCholeskyFactorSymbolic:Fill negative"); 1685*bbf0e169SBarry Smith if (!fact) SETERRQ(1,"MatIncompleteCholeskyFactorSymbolic:Missing fact argument"); 1686*bbf0e169SBarry Smith if (!mat->ops.incompletecholeskyfactorsymbolic) 1687*bbf0e169SBarry Smith SETERRQ(PETSC_ERR_SUP,"MatIncompleteCholeskyFactorSymbolic"); 1688*bbf0e169SBarry Smith if (!mat->assembled) 1689*bbf0e169SBarry Smith SETERRQ(1,"MatIncompleteCholeskyFactorSymbolic:Not for unassembled matrix"); 1690*bbf0e169SBarry Smith 1691*bbf0e169SBarry Smith PLogEventBegin(MAT_IncompleteCholeskyFactorSymbolic,mat,perm,0,0); 1692*bbf0e169SBarry Smith ierr = (*mat->ops.incompletecholeskyfactorsymbolic)(mat,perm,f,fill,fact);CHKERRQ(ierr); 1693*bbf0e169SBarry Smith PLogEventEnd(MAT_IncompleteCholeskyFactorSymbolic,mat,perm,0,0); 1694*bbf0e169SBarry Smith return 0; 1695*bbf0e169SBarry Smith } 1696*bbf0e169SBarry Smith 1697*bbf0e169SBarry Smith /*@C 1698*bbf0e169SBarry Smith MatGetArray - Returns a pointer to the element values in the matrix. 1699*bbf0e169SBarry Smith This routine is implementation dependent, and may not even work for 1700*bbf0e169SBarry Smith certain matrix types. You MUST call MatRestoreArray() when you no 1701*bbf0e169SBarry Smith longer need to access the array. 1702*bbf0e169SBarry Smith 1703*bbf0e169SBarry Smith Input Parameter: 1704*bbf0e169SBarry Smith . mat - the matrix 1705*bbf0e169SBarry Smith 1706*bbf0e169SBarry Smith Output Parameter: 1707*bbf0e169SBarry Smith . v - the location of the values 1708*bbf0e169SBarry Smith 1709*bbf0e169SBarry Smith Fortran Note: 1710*bbf0e169SBarry Smith The Fortran interface is slightly different from that given below. 1711*bbf0e169SBarry Smith See the Fortran chapter of the users manual and 1712*bbf0e169SBarry Smith petsc/src/mat/examples for details. 1713*bbf0e169SBarry Smith 1714*bbf0e169SBarry Smith .keywords: matrix, array, elements, values 1715*bbf0e169SBarry Smith 1716*bbf0e169SBarry Smith .seeaols: MatRestoreArray() 1717*bbf0e169SBarry Smith @*/ 1718*bbf0e169SBarry Smith int MatGetArray(Mat mat,Scalar **v) 1719*bbf0e169SBarry Smith { 1720*bbf0e169SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 1721*bbf0e169SBarry Smith if (!v) SETERRQ(1,"MatGetArray:Bad input, array pointer location"); 1722*bbf0e169SBarry Smith if (!mat->ops.getarray) SETERRQ(PETSC_ERR_SUP,"MatGetArray"); 1723*bbf0e169SBarry Smith return (*mat->ops.getarray)(mat,v); 1724*bbf0e169SBarry Smith } 1725*bbf0e169SBarry Smith 1726*bbf0e169SBarry Smith /*@C 1727*bbf0e169SBarry Smith MatRestoreArray - Restores the matrix after MatGetArray has been called. 1728*bbf0e169SBarry Smith 1729*bbf0e169SBarry Smith Input Parameter: 1730*bbf0e169SBarry Smith . mat - the matrix 1731*bbf0e169SBarry Smith . v - the location of the values 1732*bbf0e169SBarry Smith 1733*bbf0e169SBarry Smith Fortran Note: 1734*bbf0e169SBarry Smith The Fortran interface is slightly different from that given below. 1735*bbf0e169SBarry Smith See the users manual and petsc/src/mat/examples for details. 1736*bbf0e169SBarry Smith 1737*bbf0e169SBarry Smith .keywords: matrix, array, elements, values, resrore 1738*bbf0e169SBarry Smith 1739*bbf0e169SBarry Smith .seealso: MatGetArray() 1740*bbf0e169SBarry Smith @*/ 1741*bbf0e169SBarry Smith int MatRestoreArray(Mat mat,Scalar **v) 1742*bbf0e169SBarry Smith { 1743*bbf0e169SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 1744*bbf0e169SBarry Smith if (!v) SETERRQ(1,"MatRestoreArray:Bad input, array pointer location"); 1745*bbf0e169SBarry Smith if (!mat->ops.restorearray) SETERRQ(PETSC_ERR_SUP,"MatResroreArray"); 1746*bbf0e169SBarry Smith return (*mat->ops.restorearray)(mat,v); 1747*bbf0e169SBarry Smith } 1748*bbf0e169SBarry Smith 1749*bbf0e169SBarry Smith /*@C 1750*bbf0e169SBarry Smith MatGetSubMatrices - Extracts several submatrices from a matrix. If submat 1751*bbf0e169SBarry Smith points to an array of valid matrices, it may be reused. 1752*bbf0e169SBarry Smith 1753*bbf0e169SBarry Smith Input Parameters: 1754*bbf0e169SBarry Smith . mat - the matrix 1755*bbf0e169SBarry Smith . irow, icol - index sets of rows and columns to extract 1756*bbf0e169SBarry Smith . scall - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX 1757*bbf0e169SBarry Smith 1758*bbf0e169SBarry Smith Output Parameter: 1759*bbf0e169SBarry Smith . submat - the array of submatrices 1760*bbf0e169SBarry Smith 1761*bbf0e169SBarry Smith Notes: 1762*bbf0e169SBarry Smith When finished using the submatrices, the user should destroy 1763*bbf0e169SBarry Smith them with MatDestroySubMatrices(). 1764*bbf0e169SBarry Smith 1765*bbf0e169SBarry Smith .keywords: matrix, get, submatrix, submatrices 1766*bbf0e169SBarry Smith 1767*bbf0e169SBarry Smith .seealso: MatDestroyMatrices() 1768*bbf0e169SBarry Smith @*/ 1769*bbf0e169SBarry Smith int MatGetSubMatrices(Mat mat,int n,IS *irow,IS *icol,MatGetSubMatrixCall scall, 1770*bbf0e169SBarry Smith Mat **submat) 1771*bbf0e169SBarry Smith { 1772*bbf0e169SBarry Smith int ierr; 1773*bbf0e169SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 1774*bbf0e169SBarry Smith if (!mat->ops.getsubmatrices) SETERRQ(PETSC_ERR_SUP,"MatGetSubMatrices"); 1775*bbf0e169SBarry Smith if (!mat->assembled) SETERRQ(1,"MatGetSubMatrices:Not for unassembled matrix"); 1776*bbf0e169SBarry Smith 1777*bbf0e169SBarry Smith PLogEventBegin(MAT_GetSubMatrices,mat,0,0,0); 1778*bbf0e169SBarry Smith ierr = (*mat->ops.getsubmatrices)(mat,n,irow,icol,scall,submat); CHKERRQ(ierr); 1779*bbf0e169SBarry Smith PLogEventEnd(MAT_GetSubMatrices,mat,0,0,0); 1780*bbf0e169SBarry Smith 1781*bbf0e169SBarry Smith return 0; 1782*bbf0e169SBarry Smith } 1783*bbf0e169SBarry Smith 1784*bbf0e169SBarry Smith /*@C 1785*bbf0e169SBarry Smith MatDestroyMatrices - Destroys a set of matrices obtained with MatGetSubMatrices(). 1786*bbf0e169SBarry Smith 1787*bbf0e169SBarry Smith Input Parameters: 1788*bbf0e169SBarry Smith . n - the number of local matrices 1789*bbf0e169SBarry Smith . mat - the matrices 1790*bbf0e169SBarry Smith 1791*bbf0e169SBarry Smith .keywords: matrix, destroy, submatrix, submatrices 1792*bbf0e169SBarry Smith 1793*bbf0e169SBarry Smith .seealso: MatGetSubMatrices() 1794*bbf0e169SBarry Smith @*/ 1795*bbf0e169SBarry Smith int MatDestroyMatrices(int n,Mat **mat) 1796*bbf0e169SBarry Smith { 1797*bbf0e169SBarry Smith int ierr,i; 1798*bbf0e169SBarry Smith 1799*bbf0e169SBarry Smith PetscValidPointer(mat); 1800*bbf0e169SBarry Smith for ( i=0; i<n; i++ ) { 1801*bbf0e169SBarry Smith ierr = MatDestroy((*mat)[i]); CHKERRQ(ierr); 1802*bbf0e169SBarry Smith } 1803*bbf0e169SBarry Smith PetscFree(*mat); 1804*bbf0e169SBarry Smith return 0; 1805*bbf0e169SBarry Smith } 1806*bbf0e169SBarry Smith 1807*bbf0e169SBarry Smith /*@ 1808*bbf0e169SBarry Smith MatIncreaseOverlap - Given a set of submatrices indicated by index sets, 1809*bbf0e169SBarry Smith replaces the index by larger ones that represent submatrices with more 1810*bbf0e169SBarry Smith overlap. 1811*bbf0e169SBarry Smith 1812*bbf0e169SBarry Smith Input Parameters: 1813*bbf0e169SBarry Smith . mat - the matrix 1814*bbf0e169SBarry Smith . n - the number of index sets 1815*bbf0e169SBarry Smith . is - the array of pointers to index sets 1816*bbf0e169SBarry Smith . ov - the additional overlap requested 1817*bbf0e169SBarry Smith 1818*bbf0e169SBarry Smith .keywords: matrix, overlap, Schwarz 1819*bbf0e169SBarry Smith 1820*bbf0e169SBarry Smith .seealso: MatGetSubMatrices() 1821*bbf0e169SBarry Smith @*/ 1822*bbf0e169SBarry Smith int MatIncreaseOverlap(Mat mat,int n, IS *is,int ov) 1823*bbf0e169SBarry Smith { 1824*bbf0e169SBarry Smith int ierr; 1825*bbf0e169SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 1826*bbf0e169SBarry Smith if (!mat->assembled) SETERRQ(1,"MatIncreaseOverlap:Not for unassembled matrix"); 1827*bbf0e169SBarry Smith 1828*bbf0e169SBarry Smith if (ov == 0) return 0; 1829*bbf0e169SBarry Smith if (!mat->ops.increaseoverlap) SETERRQ(PETSC_ERR_SUP,"MatIncreaseOverlap"); 1830*bbf0e169SBarry Smith PLogEventBegin(MAT_IncreaseOverlap,mat,0,0,0); 1831*bbf0e169SBarry Smith ierr = (*mat->ops.increaseoverlap)(mat,n,is,ov); CHKERRQ(ierr); 1832*bbf0e169SBarry Smith PLogEventEnd(MAT_IncreaseOverlap,mat,0,0,0); 1833*bbf0e169SBarry Smith return 0; 1834*bbf0e169SBarry Smith } 1835*bbf0e169SBarry Smith 1836*bbf0e169SBarry Smith /*@ 1837*bbf0e169SBarry Smith MatPrintHelp - Prints all the options for the matrix. 1838*bbf0e169SBarry Smith 1839*bbf0e169SBarry Smith Input Parameter: 1840*bbf0e169SBarry Smith . mat - the matrix 1841*bbf0e169SBarry Smith 1842*bbf0e169SBarry Smith Options Database Keys: 1843*bbf0e169SBarry Smith $ -help, -h 1844*bbf0e169SBarry Smith 1845*bbf0e169SBarry Smith .keywords: mat, help 1846*bbf0e169SBarry Smith 1847*bbf0e169SBarry Smith .seealso: MatCreate(), MatCreateXXX() 1848*bbf0e169SBarry Smith @*/ 1849*bbf0e169SBarry Smith int MatPrintHelp(Mat mat) 1850*bbf0e169SBarry Smith { 1851*bbf0e169SBarry Smith static int called = 0; 1852*bbf0e169SBarry Smith MPI_Comm comm = mat->comm; 1853*bbf0e169SBarry Smith 1854*bbf0e169SBarry Smith if (!called) { 1855*bbf0e169SBarry Smith PetscPrintf(comm,"General matrix options:\n"); 1856*bbf0e169SBarry Smith PetscPrintf(comm," -mat_view_info : view basic matrix info during MatAssemblyEnd()\n"); 1857*bbf0e169SBarry Smith PetscPrintf(comm," -mat_view_info_detailed : view detailed matrix info during MatAssemblyEnd()\n"); 1858*bbf0e169SBarry Smith PetscPrintf(comm," -mat_view_draw : draw nonzero matrix structure during MatAssemblyEnd()\n"); 1859*bbf0e169SBarry Smith PetscPrintf(comm," -draw_pause <sec> : set seconds of display pause\n"); 1860*bbf0e169SBarry Smith PetscPrintf(comm," -display <name> : set alternate display\n"); 1861*bbf0e169SBarry Smith called = 1; 1862*bbf0e169SBarry Smith } 1863*bbf0e169SBarry Smith if (mat->ops.printhelp) (*mat->ops.printhelp)(mat); 1864*bbf0e169SBarry Smith return 0; 1865*bbf0e169SBarry Smith } 1866*bbf0e169SBarry Smith 1867*bbf0e169SBarry Smith /*@ 1868*bbf0e169SBarry Smith MatGetBlockSize - Returns the matrix block size; useful especially for the 1869*bbf0e169SBarry Smith block row and block diagonal formats. 1870*bbf0e169SBarry Smith 1871*bbf0e169SBarry Smith Input Parameter: 1872*bbf0e169SBarry Smith . mat - the matrix 1873*bbf0e169SBarry Smith 1874*bbf0e169SBarry Smith Output Parameter: 1875*bbf0e169SBarry Smith . bs - block size 1876*bbf0e169SBarry Smith 1877*bbf0e169SBarry Smith Notes: 1878*bbf0e169SBarry Smith $ block diagonal formats: MATSEQBDIAG, MATMPIBDIAG 1879*bbf0e169SBarry Smith $ block row formats: MATSEQBAIJ, MATMPIBAIJ 1880*bbf0e169SBarry Smith 1881*bbf0e169SBarry Smith .keywords: matrix, get, block, size 1882*bbf0e169SBarry Smith 1883*bbf0e169SBarry Smith .seealso: MatCreateSeqBAIJ(), MatCreateMPIBAIJ(), MatCreateSeqBDiag(), MatCreateMPIBDiag() 1884*bbf0e169SBarry Smith @*/ 1885*bbf0e169SBarry Smith int MatGetBlockSize(Mat mat,int *bs) 1886*bbf0e169SBarry Smith { 1887*bbf0e169SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 1888*bbf0e169SBarry Smith PetscValidIntPointer(bs); 1889*bbf0e169SBarry Smith if (!mat->ops.getblocksize) SETERRQ(PETSC_ERR_SUP,"MatGetBlockSize"); 1890*bbf0e169SBarry Smith return (*mat->ops.getblocksize)(mat,bs); 1891*bbf0e169SBarry Smith } 1892*bbf0e169SBarry Smith 1893*bbf0e169SBarry Smith /*@C 1894*bbf0e169SBarry Smith MatGetRowIJ - Returns the compress row storage i and j indices for sequential matrices. 1895*bbf0e169SBarry Smith EXPERTS ONLY. 1896*bbf0e169SBarry Smith 1897*bbf0e169SBarry Smith Input Parameters: 1898*bbf0e169SBarry Smith . mat - the matrix 1899*bbf0e169SBarry Smith . shift - 1 or zero indicating we want the indices starting at 0 or 1 1900*bbf0e169SBarry Smith . symmetric - PETSC_TRUE or PETSC_FALSE indicating the matrix data structure should be 1901*bbf0e169SBarry Smith symmetrized 1902*bbf0e169SBarry Smith 1903*bbf0e169SBarry Smith Output Parameters: 1904*bbf0e169SBarry Smith . n - number of rows and columns in the (possibly compressed) matrix 1905*bbf0e169SBarry Smith . ia - the row indices 1906*bbf0e169SBarry Smith . ja - the column indices 1907*bbf0e169SBarry Smith . done - PETSC_TRUE or PETSC_FALSE indicated that the values have been returned 1908*bbf0e169SBarry Smith @*/ 1909*bbf0e169SBarry Smith int MatGetRowIJ(Mat mat,int shift,PetscTruth symmetric,int *n,int **ia,int** ja,PetscTruth *done) 1910*bbf0e169SBarry Smith { 1911*bbf0e169SBarry Smith int ierr; 1912*bbf0e169SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 1913*bbf0e169SBarry Smith if (ia) PetscValidIntPointer(ia); 1914*bbf0e169SBarry Smith if (ja) PetscValidIntPointer(ja); 1915*bbf0e169SBarry Smith PetscValidIntPointer(done); 1916*bbf0e169SBarry Smith 1917*bbf0e169SBarry Smith if (!mat->ops.getrowij) *done = PETSC_FALSE; 1918*bbf0e169SBarry Smith else { 1919*bbf0e169SBarry Smith *done = PETSC_TRUE; 1920*bbf0e169SBarry Smith ierr = (*mat->ops.getrowij)(mat,shift,symmetric,n,ia,ja,done); CHKERRQ(ierr); 1921*bbf0e169SBarry Smith } 1922*bbf0e169SBarry Smith return 0; 1923*bbf0e169SBarry Smith } 1924*bbf0e169SBarry Smith 1925*bbf0e169SBarry Smith /*@C 1926*bbf0e169SBarry Smith MatGetColumnIJ - Returns the compress Column storage i and j indices for sequential matrices. 1927*bbf0e169SBarry Smith EXPERTS ONLY. 1928*bbf0e169SBarry Smith 1929*bbf0e169SBarry Smith Input Parameters: 1930*bbf0e169SBarry Smith . mat - the matrix 1931*bbf0e169SBarry Smith . shift - 1 or zero indicating we want the indices starting at 0 or 1 1932*bbf0e169SBarry Smith . symmetric - PETSC_TRUE or PETSC_FALSE indicating the matrix data structure should be 1933*bbf0e169SBarry Smith symmetrized 1934*bbf0e169SBarry Smith 1935*bbf0e169SBarry Smith Output Parameters: 1936*bbf0e169SBarry Smith . n - number of Columns and columns in the (possibly compressed) matrix 1937*bbf0e169SBarry Smith . ia - the Column indices 1938*bbf0e169SBarry Smith . ja - the column indices 1939*bbf0e169SBarry Smith . done - PETSC_TRUE or PETSC_FALSE indicated that the values have been returned 1940*bbf0e169SBarry Smith @*/ 1941*bbf0e169SBarry Smith int MatGetColumnIJ(Mat mat,int shift,PetscTruth symmetric,int *n,int **ia,int** ja,PetscTruth *done) 1942*bbf0e169SBarry Smith { 1943*bbf0e169SBarry Smith int ierr; 1944*bbf0e169SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 1945*bbf0e169SBarry Smith if (ia) PetscValidIntPointer(ia); 1946*bbf0e169SBarry Smith if (ja) PetscValidIntPointer(ja); 1947*bbf0e169SBarry Smith PetscValidIntPointer(done); 1948*bbf0e169SBarry Smith 1949*bbf0e169SBarry Smith if (!mat->ops.getcolumnij) *done = PETSC_FALSE; 1950*bbf0e169SBarry Smith else { 1951*bbf0e169SBarry Smith *done = PETSC_TRUE; 1952*bbf0e169SBarry Smith ierr = (*mat->ops.getcolumnij)(mat,shift,symmetric,n,ia,ja,done); CHKERRQ(ierr); 1953*bbf0e169SBarry Smith } 1954*bbf0e169SBarry Smith return 0; 1955*bbf0e169SBarry Smith } 1956*bbf0e169SBarry Smith 1957*bbf0e169SBarry Smith /*@C 1958*bbf0e169SBarry Smith MatRestoreRowIJ - Call after you are completed with the ia,ja indices obtained with 1959*bbf0e169SBarry Smith MatGetRowIJ(). EXPERTS ONLY. 1960*bbf0e169SBarry Smith 1961*bbf0e169SBarry Smith Input Parameters: 1962*bbf0e169SBarry Smith . mat - the matrix 1963*bbf0e169SBarry Smith . shift - 1 or zero indicating we want the indices starting at 0 or 1 1964*bbf0e169SBarry Smith . symmetric - PETSC_TRUE or PETSC_FALSE indicating the matrix data structure should be 1965*bbf0e169SBarry Smith symmetrized 1966*bbf0e169SBarry Smith 1967*bbf0e169SBarry Smith Output Parameters: 1968*bbf0e169SBarry Smith . n - size of (possibly compressed) matrix 1969*bbf0e169SBarry Smith . ia - the row indices 1970*bbf0e169SBarry Smith . ja - the column indices 1971*bbf0e169SBarry Smith . done - PETSC_TRUE or PETSC_FALSE indicated that the values have been returned 1972*bbf0e169SBarry Smith @*/ 1973*bbf0e169SBarry Smith int MatRestoreRowIJ(Mat mat,int shift,PetscTruth symmetric,int *n,int **ia,int** ja,PetscTruth *done) 1974*bbf0e169SBarry Smith { 1975*bbf0e169SBarry Smith int ierr; 1976*bbf0e169SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 1977*bbf0e169SBarry Smith if (ia) PetscValidIntPointer(ia); 1978*bbf0e169SBarry Smith if (ja) PetscValidIntPointer(ja); 1979*bbf0e169SBarry Smith PetscValidIntPointer(done); 1980*bbf0e169SBarry Smith 1981*bbf0e169SBarry Smith if (!mat->ops.restorerowij) *done = PETSC_FALSE; 1982*bbf0e169SBarry Smith else { 1983*bbf0e169SBarry Smith *done = PETSC_TRUE; 1984*bbf0e169SBarry Smith ierr = (*mat->ops.restorerowij)(mat,shift,symmetric,n,ia,ja,done); CHKERRQ(ierr); 1985*bbf0e169SBarry Smith } 1986*bbf0e169SBarry Smith return 0; 1987*bbf0e169SBarry Smith } 1988*bbf0e169SBarry Smith 1989*bbf0e169SBarry Smith /*@C 1990*bbf0e169SBarry Smith MatRestoreColumnIJ - Call after you are completed with the ia,ja indices obtained with 1991*bbf0e169SBarry Smith MatGetColumnIJ(). EXPERTS ONLY. 1992*bbf0e169SBarry Smith 1993*bbf0e169SBarry Smith Input Parameters: 1994*bbf0e169SBarry Smith . mat - the matrix 1995*bbf0e169SBarry Smith . shift - 1 or zero indicating we want the indices starting at 0 or 1 1996*bbf0e169SBarry Smith . symmetric - PETSC_TRUE or PETSC_FALSE indicating the matrix data structure should be 1997*bbf0e169SBarry Smith symmetrized 1998*bbf0e169SBarry Smith 1999*bbf0e169SBarry Smith Output Parameters: 2000*bbf0e169SBarry Smith . n - size of (possibly compressed) matrix 2001*bbf0e169SBarry Smith . ia - the Column indices 2002*bbf0e169SBarry Smith . ja - the column indices 2003*bbf0e169SBarry Smith . done - PETSC_TRUE or PETSC_FALSE indicated that the values have been returned 2004*bbf0e169SBarry Smith @*/ 2005*bbf0e169SBarry Smith int MatRestoreColumnIJ(Mat mat,int shift,PetscTruth symmetric,int *n,int **ia,int** ja,PetscTruth *done) 2006*bbf0e169SBarry Smith { 2007*bbf0e169SBarry Smith int ierr; 2008*bbf0e169SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 2009*bbf0e169SBarry Smith if (ia) PetscValidIntPointer(ia); 2010*bbf0e169SBarry Smith if (ja) PetscValidIntPointer(ja); 2011*bbf0e169SBarry Smith PetscValidIntPointer(done); 2012*bbf0e169SBarry Smith 2013*bbf0e169SBarry Smith if (!mat->ops.restorecolumnij) *done = PETSC_FALSE; 2014*bbf0e169SBarry Smith else { 2015*bbf0e169SBarry Smith *done = PETSC_TRUE; 2016*bbf0e169SBarry Smith ierr = (*mat->ops.restorecolumnij)(mat,shift,symmetric,n,ia,ja,done); CHKERRQ(ierr); 2017*bbf0e169SBarry Smith } 2018*bbf0e169SBarry Smith return 0; 2019*bbf0e169SBarry Smith } 2020*bbf0e169SBarry Smith 2021*bbf0e169SBarry Smith 2022