#define PETSCMAT_DLL /* This is where the abstract matrix operations are defined */ #include "include/private/matimpl.h" /*I "petscmat.h" I*/ #include "private/vecimpl.h" /* Logging support */ PetscCookie PETSCMAT_DLLEXPORT MAT_COOKIE = 0; PetscEvent MAT_Mult = 0, MAT_Mults = 0, MAT_MultConstrained = 0, MAT_MultAdd = 0, MAT_MultTranspose = 0; PetscEvent MAT_MultTransposeConstrained = 0, MAT_MultTransposeAdd = 0, MAT_Solve = 0, MAT_Solves = 0, MAT_SolveAdd = 0, MAT_SolveTranspose = 0, MAT_MatSolve = 0; PetscEvent MAT_SolveTransposeAdd = 0, MAT_Relax = 0, MAT_ForwardSolve = 0, MAT_BackwardSolve = 0, MAT_LUFactor = 0, MAT_LUFactorSymbolic = 0; PetscEvent MAT_LUFactorNumeric = 0, MAT_CholeskyFactor = 0, MAT_CholeskyFactorSymbolic = 0, MAT_CholeskyFactorNumeric = 0, MAT_ILUFactor = 0; PetscEvent MAT_ILUFactorSymbolic = 0, MAT_ICCFactorSymbolic = 0, MAT_Copy = 0, MAT_Convert = 0, MAT_Scale = 0, MAT_AssemblyBegin = 0; PetscEvent MAT_AssemblyEnd = 0, MAT_SetValues = 0, MAT_GetValues = 0, MAT_GetRow = 0, MAT_GetRowIJ = 0, MAT_GetSubMatrices = 0, MAT_GetColoring = 0, MAT_GetOrdering = 0, MAT_GetRedundantMatrix = 0; PetscEvent MAT_IncreaseOverlap = 0, MAT_Partitioning = 0, MAT_ZeroEntries = 0, MAT_Load = 0, MAT_View = 0, MAT_AXPY = 0, MAT_FDColoringCreate = 0; PetscEvent MAT_FDColoringApply = 0,MAT_Transpose = 0,MAT_FDColoringFunction = 0; PetscEvent MAT_MatMult = 0, MAT_MatMultSymbolic = 0, MAT_MatMultNumeric = 0; PetscEvent MAT_PtAP = 0, MAT_PtAPSymbolic = 0, MAT_PtAPNumeric = 0; PetscEvent MAT_MatMultTranspose = 0, MAT_MatMultTransposeSymbolic = 0, MAT_MatMultTransposeNumeric = 0; /* nasty global values for MatSetValue() */ PetscInt PETSCMAT_DLLEXPORT MatSetValue_Row = 0; PetscInt PETSCMAT_DLLEXPORT MatSetValue_Column = 0; PetscScalar PETSCMAT_DLLEXPORT MatSetValue_Value = 0.0; #undef __FUNCT__ #define __FUNCT__ "MatRealPart" /*@ MatRealPart - Zeros out the imaginary part of the matrix Collective on Mat Input Parameters: . mat - the matrix Level: advanced .seealso: MatImaginaryPart() @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatRealPart(Mat mat) { PetscErrorCode ierr; PetscFunctionBegin; PetscValidHeaderSpecific(mat,MAT_COOKIE,1); PetscValidType(mat,1); if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); if (!mat->ops->realpart) SETERRQ1(PETSC_ERR_SUP,"Mat type %s",mat->type_name); ierr = MatPreallocated(mat);CHKERRQ(ierr); ierr = (*mat->ops->realpart)(mat);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatImaginaryPart" /*@ MatImaginaryPart - Moves the imaginary part of the matrix to the real part and zeros the imaginary part Collective on Mat Input Parameters: . mat - the matrix Level: advanced .seealso: MatRealPart() @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatImaginaryPart(Mat mat) { PetscErrorCode ierr; PetscFunctionBegin; PetscValidHeaderSpecific(mat,MAT_COOKIE,1); PetscValidType(mat,1); if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); if (!mat->ops->imaginarypart) SETERRQ1(PETSC_ERR_SUP,"Mat type %s",mat->type_name); ierr = MatPreallocated(mat);CHKERRQ(ierr); ierr = (*mat->ops->imaginarypart)(mat);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatGetRow" /*@C MatGetRow - Gets a row of a matrix. You MUST call MatRestoreRow() for each row that you get to ensure that your application does not bleed memory. Not Collective Input Parameters: + mat - the matrix - row - the row to get Output Parameters: + ncols - if not NULL, the number of nonzeros in the row . cols - if not NULL, the column numbers - vals - if not NULL, the values Notes: This routine is provided for people who need to have direct access to the structure of a matrix. We hope that we provide enough high-level matrix routines that few users will need it. MatGetRow() always returns 0-based column indices, regardless of whether the internal representation is 0-based (default) or 1-based. For better efficiency, set cols and/or vals to PETSC_NULL if you do not wish to extract these quantities. The user can only examine the values extracted with MatGetRow(); the values cannot be altered. To change the matrix entries, one must use MatSetValues(). You can only have one call to MatGetRow() outstanding for a particular matrix at a time, per processor. MatGetRow() can only obtain rows associated with the given processor, it cannot get rows from the other processors; for that we suggest using MatGetSubMatrices(), then MatGetRow() on the submatrix. The row indix passed to MatGetRows() is in the global number of rows. Fortran Notes: The calling sequence from Fortran is .vb MatGetRow(matrix,row,ncols,cols,values,ierr) Mat matrix (input) integer row (input) integer ncols (output) integer cols(maxcols) (output) double precision (or double complex) values(maxcols) output .ve where maxcols >= maximum nonzeros in any row of the matrix. Caution: Do not try to change the contents of the output arrays (cols and vals). In some cases, this may corrupt the matrix. Level: advanced Concepts: matrices^row access .seealso: MatRestoreRow(), MatSetValues(), MatGetValues(), MatGetSubMatrices(), MatGetDiagonal() @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatGetRow(Mat mat,PetscInt row,PetscInt *ncols,const PetscInt *cols[],const PetscScalar *vals[]) { PetscErrorCode ierr; PetscInt incols; PetscFunctionBegin; PetscValidHeaderSpecific(mat,MAT_COOKIE,1); PetscValidType(mat,1); if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); if (!mat->ops->getrow) SETERRQ1(PETSC_ERR_SUP,"Mat type %s",mat->type_name); ierr = MatPreallocated(mat);CHKERRQ(ierr); ierr = PetscLogEventBegin(MAT_GetRow,mat,0,0,0);CHKERRQ(ierr); ierr = (*mat->ops->getrow)(mat,row,&incols,(PetscInt **)cols,(PetscScalar **)vals);CHKERRQ(ierr); if (ncols) *ncols = incols; ierr = PetscLogEventEnd(MAT_GetRow,mat,0,0,0);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatConjugate" /*@ MatConjugate - replaces the matrix values with their complex conjugates Collective on Mat Input Parameters: . mat - the matrix Level: advanced .seealso: VecConjugate() @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatConjugate(Mat mat) { PetscErrorCode ierr; PetscFunctionBegin; PetscValidHeaderSpecific(mat,MAT_COOKIE,1); if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); if (!mat->ops->conjugate) SETERRQ(PETSC_ERR_SUP,"Not provided for this matrix format, send email to petsc-maint@mcs.anl.gov"); ierr = (*mat->ops->conjugate)(mat);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatRestoreRow" /*@C MatRestoreRow - Frees any temporary space allocated by MatGetRow(). Not Collective Input Parameters: + mat - the matrix . row - the row to get . ncols, cols - the number of nonzeros and their columns - vals - if nonzero the column values Notes: This routine should be called after you have finished examining the entries. Fortran Notes: The calling sequence from Fortran is .vb MatRestoreRow(matrix,row,ncols,cols,values,ierr) Mat matrix (input) integer row (input) integer ncols (output) integer cols(maxcols) (output) double precision (or double complex) values(maxcols) output .ve Where maxcols >= maximum nonzeros in any row of the matrix. In Fortran MatRestoreRow() MUST be called after MatGetRow() before another call to MatGetRow() can be made. Level: advanced .seealso: MatGetRow() @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRow(Mat mat,PetscInt row,PetscInt *ncols,const PetscInt *cols[],const PetscScalar *vals[]) { PetscErrorCode ierr; PetscFunctionBegin; PetscValidHeaderSpecific(mat,MAT_COOKIE,1); PetscValidIntPointer(ncols,3); if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); if (!mat->ops->restorerow) PetscFunctionReturn(0); ierr = (*mat->ops->restorerow)(mat,row,ncols,(PetscInt **)cols,(PetscScalar **)vals);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatGetRowUpperTriangular" /*@C MatGetRowUpperTriangular - Sets a flag to enable calls to MatGetRow() for matrix in MATSBAIJ format. You should call MatRestoreRowUpperTriangular() after calling MatGetRow/MatRestoreRow() to disable the flag. Not Collective Input Parameters: + mat - the matrix Notes: The flag is to ensure that users are aware of MatGetRow() only provides the upper trianglular part of the row for the matrices in MATSBAIJ format. Level: advanced Concepts: matrices^row access .seealso: MatRestoreRowRowUpperTriangular() @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowUpperTriangular(Mat mat) { PetscErrorCode ierr; PetscFunctionBegin; PetscValidHeaderSpecific(mat,MAT_COOKIE,1); PetscValidType(mat,1); if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); if (!mat->ops->getrowuppertriangular) SETERRQ1(PETSC_ERR_SUP,"Mat type %s",mat->type_name); ierr = MatPreallocated(mat);CHKERRQ(ierr); ierr = (*mat->ops->getrowuppertriangular)(mat);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatRestoreRowUpperTriangular" /*@C MatRestoreRowUpperTriangular - Disable calls to MatGetRow() for matrix in MATSBAIJ format. Not Collective Input Parameters: + mat - the matrix Notes: This routine should be called after you have finished MatGetRow/MatRestoreRow(). Level: advanced .seealso: MatGetRowUpperTriangular() @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRowUpperTriangular(Mat mat) { PetscErrorCode ierr; PetscFunctionBegin; PetscValidHeaderSpecific(mat,MAT_COOKIE,1); if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); if (!mat->ops->restorerowuppertriangular) PetscFunctionReturn(0); ierr = (*mat->ops->restorerowuppertriangular)(mat);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatSetOptionsPrefix" /*@C MatSetOptionsPrefix - Sets the prefix used for searching for all Mat options in the database. Collective on Mat Input Parameter: + A - the Mat context - prefix - the prefix to prepend to all option names Notes: A hyphen (-) must NOT be given at the beginning of the prefix name. The first character of all runtime options is AUTOMATICALLY the hyphen. Level: advanced .keywords: Mat, set, options, prefix, database .seealso: MatSetFromOptions() @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatSetOptionsPrefix(Mat A,const char prefix[]) { PetscErrorCode ierr; PetscFunctionBegin; PetscValidHeaderSpecific(A,MAT_COOKIE,1); ierr = PetscObjectSetOptionsPrefix((PetscObject)A,prefix);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatAppendOptionsPrefix" /*@C MatAppendOptionsPrefix - Appends to the prefix used for searching for all Mat options in the database. Collective on Mat Input Parameters: + A - the Mat context - prefix - the prefix to prepend to all option names Notes: A hyphen (-) must NOT be given at the beginning of the prefix name. The first character of all runtime options is AUTOMATICALLY the hyphen. Level: advanced .keywords: Mat, append, options, prefix, database .seealso: MatGetOptionsPrefix() @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatAppendOptionsPrefix(Mat A,const char prefix[]) { PetscErrorCode ierr; PetscFunctionBegin; PetscValidHeaderSpecific(A,MAT_COOKIE,1); ierr = PetscObjectAppendOptionsPrefix((PetscObject)A,prefix);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatGetOptionsPrefix" /*@C MatGetOptionsPrefix - Sets the prefix used for searching for all Mat options in the database. Not Collective Input Parameter: . A - the Mat context Output Parameter: . prefix - pointer to the prefix string used Notes: On the fortran side, the user should pass in a string 'prefix' of sufficient length to hold the prefix. Level: advanced .keywords: Mat, get, options, prefix, database .seealso: MatAppendOptionsPrefix() @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatGetOptionsPrefix(Mat A,const char *prefix[]) { PetscErrorCode ierr; PetscFunctionBegin; PetscValidHeaderSpecific(A,MAT_COOKIE,1); ierr = PetscObjectGetOptionsPrefix((PetscObject)A,prefix);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatSetUp" /*@ MatSetUp - Sets up the internal matrix data structures for the later use. Collective on Mat Input Parameters: . A - the Mat context Notes: For basic use of the Mat classes the user need not explicitly call MatSetUp(), since these actions will happen automatically. Level: advanced .keywords: Mat, setup .seealso: MatCreate(), MatDestroy() @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatSetUp(Mat A) { PetscErrorCode ierr; PetscFunctionBegin; PetscValidHeaderSpecific(A,MAT_COOKIE,1); ierr = MatSetUpPreallocation(A);CHKERRQ(ierr); ierr = MatSetFromOptions(A);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatView" /*@C MatView - Visualizes a matrix object. Collective on Mat Input Parameters: + mat - the matrix - viewer - visualization context Notes: The available visualization contexts include + PETSC_VIEWER_STDOUT_SELF - standard output (default) . PETSC_VIEWER_STDOUT_WORLD - synchronized standard output where only the first processor opens the file. All other processors send their data to the first processor to print. - PETSC_VIEWER_DRAW_WORLD - graphical display of nonzero structure The user can open alternative visualization contexts with + PetscViewerASCIIOpen() - Outputs matrix to a specified file . PetscViewerBinaryOpen() - Outputs matrix in binary to a specified file; corresponding input uses MatLoad() . PetscViewerDrawOpen() - Outputs nonzero matrix structure to an X window display - PetscViewerSocketOpen() - Outputs matrix to Socket viewer. Currently only the sequential dense and AIJ matrix types support the Socket viewer. The user can call PetscViewerSetFormat() to specify the output format of ASCII printed objects (when using PETSC_VIEWER_STDOUT_SELF, PETSC_VIEWER_STDOUT_WORLD and PetscViewerASCIIOpen). Available formats include + PETSC_VIEWER_ASCII_DEFAULT - default, prints matrix contents . PETSC_VIEWER_ASCII_MATLAB - prints matrix contents in Matlab format . PETSC_VIEWER_ASCII_DENSE - prints entire matrix including zeros . PETSC_VIEWER_ASCII_COMMON - prints matrix contents, using a sparse format common among all matrix types . PETSC_VIEWER_ASCII_IMPL - prints matrix contents, using an implementation-specific format (which is in many cases the same as the default) . PETSC_VIEWER_ASCII_INFO - prints basic information about the matrix size and structure (not the matrix entries) . PETSC_VIEWER_ASCII_INFO_DETAIL - prints more detailed information about the matrix structure Options Database Keys: + -mat_view_info - Prints info on matrix at conclusion of MatEndAssembly() . -mat_view_info_detailed - Prints more detailed info . -mat_view - Prints matrix in ASCII format . -mat_view_matlab - Prints matrix in Matlab format . -mat_view_draw - PetscDraws nonzero structure of matrix, using MatView() and PetscDrawOpenX(). . -display - Sets display name (default is host) . -draw_pause - Sets number of seconds to pause after display . -mat_view_socket - Sends matrix to socket, can be accessed from Matlab (see users manual) . -viewer_socket_machine . -viewer_socket_port . -mat_view_binary - save matrix to file in binary format - -viewer_binary_filename Level: beginner Concepts: matrices^viewing Concepts: matrices^plotting Concepts: matrices^printing .seealso: PetscViewerSetFormat(), PetscViewerASCIIOpen(), PetscViewerDrawOpen(), PetscViewerSocketOpen(), PetscViewerBinaryOpen(), MatLoad() @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatView(Mat mat,PetscViewer viewer) { PetscErrorCode ierr; PetscInt rows,cols; PetscTruth iascii; const char *cstr; PetscViewerFormat format; PetscFunctionBegin; PetscValidHeaderSpecific(mat,MAT_COOKIE,1); PetscValidType(mat,1); if (!viewer) { ierr = PetscViewerASCIIGetStdout(mat->comm,&viewer);CHKERRQ(ierr); } PetscValidHeaderSpecific(viewer,PETSC_VIEWER_COOKIE,2); PetscCheckSameComm(mat,1,viewer,2); if (!mat->assembled) SETERRQ(PETSC_ERR_ORDER,"Must call MatAssemblyBegin/End() before viewing matrix"); ierr = MatPreallocated(mat);CHKERRQ(ierr); ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); if (iascii) { ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) { if (mat->prefix) { ierr = PetscViewerASCIIPrintf(viewer,"Matrix Object:(%s)\n",mat->prefix);CHKERRQ(ierr); } else { ierr = PetscViewerASCIIPrintf(viewer,"Matrix Object:\n");CHKERRQ(ierr); } ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); ierr = MatGetType(mat,&cstr);CHKERRQ(ierr); ierr = MatGetSize(mat,&rows,&cols);CHKERRQ(ierr); ierr = PetscViewerASCIIPrintf(viewer,"type=%s, rows=%D, cols=%D\n",cstr,rows,cols);CHKERRQ(ierr); if (mat->ops->getinfo) { MatInfo info; ierr = MatGetInfo(mat,MAT_GLOBAL_SUM,&info);CHKERRQ(ierr); ierr = PetscViewerASCIIPrintf(viewer,"total: nonzeros=%D, allocated nonzeros=%D\n", (PetscInt)info.nz_used,(PetscInt)info.nz_allocated);CHKERRQ(ierr); } } } if (mat->ops->view) { ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); ierr = (*mat->ops->view)(mat,viewer);CHKERRQ(ierr); ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); } else if (!iascii) { SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported",((PetscObject)viewer)->type_name); } if (iascii) { ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) { ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); } } PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatScaleSystem" /*@C MatScaleSystem - Scale a vector solution and right hand side to match the scaling of a scaled matrix. Collective on Mat Input Parameter: + mat - the matrix . b - right hand side vector (or PETSC_NULL) - x - solution vector (or PETSC_NULL) Notes: For AIJ, BAIJ, and BDiag matrix formats, the matrices are not internally scaled, so this does nothing. For MPIROWBS it permutes and diagonally scales. The KSP methods automatically call this routine when required (via PCPreSolve()) so it is rarely used directly. Level: Developer Concepts: matrices^scaling .seealso: MatUseScaledForm(), MatUnScaleSystem() @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatScaleSystem(Mat mat,Vec b,Vec x) { PetscErrorCode ierr; PetscFunctionBegin; PetscValidHeaderSpecific(mat,MAT_COOKIE,1); PetscValidType(mat,1); ierr = MatPreallocated(mat);CHKERRQ(ierr); if (x) {PetscValidHeaderSpecific(x,VEC_COOKIE,2);PetscCheckSameComm(mat,1,x,2);} if (b) {PetscValidHeaderSpecific(b,VEC_COOKIE,3);PetscCheckSameComm(mat,1,b,3);} if (mat->ops->scalesystem) { ierr = (*mat->ops->scalesystem)(mat,b,x);CHKERRQ(ierr); } PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatUnScaleSystem" /*@C MatUnScaleSystem - Unscales a vector solution and right hand side to match the original scaling of a scaled matrix. Collective on Mat Input Parameter: + mat - the matrix . b - right hand side vector (or PETSC_NULL) - x - solution vector (or PETSC_NULL) Notes: For AIJ, BAIJ, and BDiag matrix formats, the matrices are not internally scaled, so this does nothing. For MPIROWBS it permutes and diagonally scales. The KSP methods automatically call this routine when required (via PCPreSolve()) so it is rarely used directly. Level: Developer .seealso: MatUseScaledForm(), MatScaleSystem() @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatUnScaleSystem(Mat mat,Vec b,Vec x) { PetscErrorCode ierr; PetscFunctionBegin; PetscValidHeaderSpecific(mat,MAT_COOKIE,1); PetscValidType(mat,1); ierr = MatPreallocated(mat);CHKERRQ(ierr); if (x) {PetscValidHeaderSpecific(x,VEC_COOKIE,2);PetscCheckSameComm(mat,1,x,2);} if (b) {PetscValidHeaderSpecific(b,VEC_COOKIE,3);PetscCheckSameComm(mat,1,b,3);} if (mat->ops->unscalesystem) { ierr = (*mat->ops->unscalesystem)(mat,b,x);CHKERRQ(ierr); } PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatUseScaledForm" /*@C MatUseScaledForm - For matrix storage formats that scale the matrix (for example MPIRowBS matrices are diagonally scaled on assembly) indicates matrix operations (MatMult() etc) are applied using the scaled matrix. Collective on Mat Input Parameter: + mat - the matrix - scaled - PETSC_TRUE for applying the scaled, PETSC_FALSE for applying the original matrix Notes: For scaled matrix formats, applying the original, unscaled matrix will be slightly more expensive Level: Developer .seealso: MatScaleSystem(), MatUnScaleSystem() @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatUseScaledForm(Mat mat,PetscTruth scaled) { PetscErrorCode ierr; PetscFunctionBegin; PetscValidHeaderSpecific(mat,MAT_COOKIE,1); PetscValidType(mat,1); ierr = MatPreallocated(mat);CHKERRQ(ierr); if (mat->ops->usescaledform) { ierr = (*mat->ops->usescaledform)(mat,scaled);CHKERRQ(ierr); } PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatDestroy" /*@ MatDestroy - Frees space taken by a matrix. Collective on Mat Input Parameter: . A - the matrix Level: beginner @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatDestroy(Mat A) { PetscErrorCode ierr; PetscFunctionBegin; PetscValidHeaderSpecific(A,MAT_COOKIE,1); if (--A->refct > 0) PetscFunctionReturn(0); ierr = MatPreallocated(A);CHKERRQ(ierr); /* if memory was published with AMS then destroy it */ ierr = PetscObjectDepublish(A);CHKERRQ(ierr); if (A->ops->destroy) { ierr = (*A->ops->destroy)(A);CHKERRQ(ierr); } if (A->mapping) { ierr = ISLocalToGlobalMappingDestroy(A->mapping);CHKERRQ(ierr); } if (A->bmapping) { ierr = ISLocalToGlobalMappingDestroy(A->bmapping);CHKERRQ(ierr); } ierr = PetscFree(A->rmap.range);CHKERRQ(ierr); ierr = PetscFree(A->cmap.range);CHKERRQ(ierr); if (A->spptr){ierr = PetscFree(A->spptr);CHKERRQ(ierr);} ierr = PetscHeaderDestroy(A);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatValid" /*@ MatValid - Checks whether a matrix object is valid. Collective on Mat Input Parameter: . m - the matrix to check Output Parameter: flg - flag indicating matrix status, either PETSC_TRUE if matrix is valid, or PETSC_FALSE otherwise. Level: developer Concepts: matrices^validity @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatValid(Mat m,PetscTruth *flg) { PetscFunctionBegin; PetscValidIntPointer(flg,1); if (!m) *flg = PETSC_FALSE; else if (m->cookie != MAT_COOKIE) *flg = PETSC_FALSE; else *flg = PETSC_TRUE; PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatSetValues" /*@ MatSetValues - Inserts or adds a block of values into a matrix. These values may be cached, so MatAssemblyBegin() and MatAssemblyEnd() MUST be called after all calls to MatSetValues() have been completed. Not Collective Input Parameters: + mat - the matrix . v - a logically two-dimensional array of values . m, idxm - the number of rows and their global indices . n, idxn - the number of columns and their global indices - addv - either ADD_VALUES or INSERT_VALUES, where ADD_VALUES adds values to any existing entries, and INSERT_VALUES replaces existing entries with new values Notes: By default the values, v, are row-oriented and unsorted. See MatSetOption() for other options. Calls to MatSetValues() with the INSERT_VALUES and ADD_VALUES options cannot be mixed without intervening calls to the assembly routines. MatSetValues() uses 0-based row and column numbers in Fortran as well as in C. Negative indices may be passed in idxm and idxn, these rows and columns are simply ignored. This allows easily inserting element stiffness matrices with homogeneous Dirchlet boundary conditions that you don't want represented in the matrix. Efficiency Alert: The routine MatSetValuesBlocked() may offer much better efficiency for users of block sparse formats (MATSEQBAIJ and MATMPIBAIJ). Level: beginner Concepts: matrices^putting entries in .seealso: MatSetOption(), MatAssemblyBegin(), MatAssemblyEnd(), MatSetValuesBlocked(), MatSetValuesLocal(), InsertMode, INSERT_VALUES, ADD_VALUES @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatSetValues(Mat mat,PetscInt m,const PetscInt idxm[],PetscInt n,const PetscInt idxn[],const PetscScalar v[],InsertMode addv) { PetscErrorCode ierr; PetscFunctionBegin; if (!m || !n) PetscFunctionReturn(0); /* no values to insert */ PetscValidHeaderSpecific(mat,MAT_COOKIE,1); PetscValidType(mat,1); PetscValidIntPointer(idxm,3); PetscValidIntPointer(idxn,5); PetscValidScalarPointer(v,6); ierr = MatPreallocated(mat);CHKERRQ(ierr); if (mat->insertmode == NOT_SET_VALUES) { mat->insertmode = addv; } #if defined(PETSC_USE_DEBUG) else if (mat->insertmode != addv) { SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Cannot mix add values and insert values"); } if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); #endif if (mat->assembled) { mat->was_assembled = PETSC_TRUE; mat->assembled = PETSC_FALSE; } ierr = PetscLogEventBegin(MAT_SetValues,mat,0,0,0);CHKERRQ(ierr); if (!mat->ops->setvalues) SETERRQ1(PETSC_ERR_SUP,"Mat type %s",mat->type_name); ierr = (*mat->ops->setvalues)(mat,m,idxm,n,idxn,v,addv);CHKERRQ(ierr); ierr = PetscLogEventEnd(MAT_SetValues,mat,0,0,0);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatSetValuesRowLocal" /*@ MatSetValuesRowLocal - Inserts a row (block row for BAIJ matrices) of nonzero values into a matrix Not Collective Input Parameters: + mat - the matrix . row - the (block) row to set - v - a logically two-dimensional array of values Notes: By the values, v, are column-oriented (for the block version) and sorted All the nonzeros in the row must be provided The matrix must have previously had its column indices set The row must belong to this process Level: intermediate Concepts: matrices^putting entries in .seealso: MatSetOption(), MatAssemblyBegin(), MatAssemblyEnd(), MatSetValuesBlocked(), MatSetValuesLocal(), InsertMode, INSERT_VALUES, ADD_VALUES, MatSetValues(), MatSetValuesRow(), MatSetLocalToGlobalMapping() @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesRowLocal(Mat mat,PetscInt row,const PetscScalar v[]) { PetscErrorCode ierr; PetscFunctionBegin; PetscValidHeaderSpecific(mat,MAT_COOKIE,1); PetscValidType(mat,1); PetscValidScalarPointer(v,2); ierr = MatSetValuesRow(mat, mat->mapping->indices[row],v);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatSetValuesRow" /*@ MatSetValuesRow - Inserts a row (block row for BAIJ matrices) of nonzero values into a matrix Not Collective Input Parameters: + mat - the matrix . row - the (block) row to set - v - a logically two-dimensional array of values Notes: By the values, v, are column-oriented (for the block version) and sorted All the nonzeros in the row must be provided The matrix must have previously had its column indices set The row must belong to this process Level: intermediate Concepts: matrices^putting entries in .seealso: MatSetOption(), MatAssemblyBegin(), MatAssemblyEnd(), MatSetValuesBlocked(), MatSetValuesLocal(), InsertMode, INSERT_VALUES, ADD_VALUES, MatSetValues() @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesRow(Mat mat,PetscInt row,const PetscScalar v[]) { PetscErrorCode ierr; PetscFunctionBegin; PetscValidHeaderSpecific(mat,MAT_COOKIE,1); PetscValidType(mat,1); PetscValidScalarPointer(v,2); #if defined(PETSC_USE_DEBUG) if (mat->insertmode == ADD_VALUES) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Cannot mix add and insert values"); if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); #endif mat->insertmode = INSERT_VALUES; if (mat->assembled) { mat->was_assembled = PETSC_TRUE; mat->assembled = PETSC_FALSE; } ierr = PetscLogEventBegin(MAT_SetValues,mat,0,0,0);CHKERRQ(ierr); if (!mat->ops->setvaluesrow) SETERRQ1(PETSC_ERR_SUP,"Mat type %s",mat->type_name); ierr = (*mat->ops->setvaluesrow)(mat,row,v);CHKERRQ(ierr); ierr = PetscLogEventEnd(MAT_SetValues,mat,0,0,0);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatSetValuesStencil" /*@ MatSetValuesStencil - Inserts or adds a block of values into a matrix. Using structured grid indexing Not Collective Input Parameters: + mat - the matrix . v - a logically two-dimensional array of values . m - number of rows being entered . idxm - grid coordinates (and component number when dof > 1) for matrix rows being entered . n - number of columns being entered . idxn - grid coordinates (and component number when dof > 1) for matrix columns being entered - addv - either ADD_VALUES or INSERT_VALUES, where ADD_VALUES adds values to any existing entries, and INSERT_VALUES replaces existing entries with new values Notes: By default the values, v, are row-oriented and unsorted. See MatSetOption() for other options. Calls to MatSetValuesStencil() with the INSERT_VALUES and ADD_VALUES options cannot be mixed without intervening calls to the assembly routines. The grid coordinates are across the entire grid, not just the local portion MatSetValuesStencil() uses 0-based row and column numbers in Fortran as well as in C. For setting/accessing vector values via array coordinates you can use the DAVecGetArray() routine In order to use this routine you must either obtain the matrix with DAGetMatrix() or call MatSetLocalToGlobalMapping() and MatSetStencil() first. The columns and rows in the stencil passed in MUST be contained within the ghost region of the given process as set with DACreateXXX() or MatSetStencil(). For example, if you create a DA with an overlap of one grid level and on a particular process its first local nonghost x logical coordinate is 6 (so its first ghost x logical coordinate is 5) the first i index you can use in your column and row indices in MatSetStencil() is 5. In Fortran idxm and idxn should be declared as $ MatStencil idxm(4,m),idxn(4,n) and the values inserted using $ idxm(MatStencil_i,1) = i $ idxm(MatStencil_j,1) = j $ idxm(MatStencil_k,1) = k $ idxm(MatStencil_c,1) = c etc For periodic boundary conditions use negative indices for values to the left (below 0; that are to be obtained by wrapping values from right edge). For values to the right of the last entry using that index plus one etc to obtain values that obtained by wrapping the values from the left edge. For indices that don't mean anything for your case (like the k index when working in 2d) or the c index when you have a single value per point) you can skip filling those indices. Inspired by the structured grid interface to the HYPRE package (http://www.llnl.gov/CASC/hypre) Efficiency Alert: The routine MatSetValuesBlockedStencil() may offer much better efficiency for users of block sparse formats (MATSEQBAIJ and MATMPIBAIJ). Level: beginner Concepts: matrices^putting entries in .seealso: MatSetOption(), MatAssemblyBegin(), MatAssemblyEnd(), MatSetValuesBlocked(), MatSetValuesLocal() MatSetValues(), MatSetValuesBlockedStencil(), MatSetStencil(), DAGetMatrix(), DAVecGetArray(), MatStencil @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesStencil(Mat mat,PetscInt m,const MatStencil idxm[],PetscInt n,const MatStencil idxn[],const PetscScalar v[],InsertMode addv) { PetscErrorCode ierr; PetscInt j,i,jdxm[128],jdxn[256],dim = mat->stencil.dim,*dims = mat->stencil.dims+1,tmp; PetscInt *starts = mat->stencil.starts,*dxm = (PetscInt*)idxm,*dxn = (PetscInt*)idxn,sdim = dim - (1 - (PetscInt)mat->stencil.noc); PetscFunctionBegin; if (!m || !n) PetscFunctionReturn(0); /* no values to insert */ PetscValidHeaderSpecific(mat,MAT_COOKIE,1); PetscValidType(mat,1); PetscValidIntPointer(idxm,3); PetscValidIntPointer(idxn,5); PetscValidScalarPointer(v,6); if (m > 128) SETERRQ1(PETSC_ERR_SUP,"Can only set 128 rows at a time; trying to set %D",m); if (n > 256) SETERRQ1(PETSC_ERR_SUP,"Can only set 256 columns at a time; trying to set %D",n); for (i=0; istencil.noc) dxm++; jdxm[i] = tmp; } for (i=0; istencil.noc) dxn++; jdxn[i] = tmp; } ierr = MatSetValuesLocal(mat,m,jdxm,n,jdxn,v,addv);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatSetValuesBlockedStencil" /*@C MatSetValuesBlockedStencil - Inserts or adds a block of values into a matrix. Using structured grid indexing Not Collective Input Parameters: + mat - the matrix . v - a logically two-dimensional array of values . m - number of rows being entered . idxm - grid coordinates for matrix rows being entered . n - number of columns being entered . idxn - grid coordinates for matrix columns being entered - addv - either ADD_VALUES or INSERT_VALUES, where ADD_VALUES adds values to any existing entries, and INSERT_VALUES replaces existing entries with new values Notes: By default the values, v, are row-oriented and unsorted. See MatSetOption() for other options. Calls to MatSetValuesBlockedStencil() with the INSERT_VALUES and ADD_VALUES options cannot be mixed without intervening calls to the assembly routines. The grid coordinates are across the entire grid, not just the local portion MatSetValuesBlockedStencil() uses 0-based row and column numbers in Fortran as well as in C. For setting/accessing vector values via array coordinates you can use the DAVecGetArray() routine In order to use this routine you must either obtain the matrix with DAGetMatrix() or call MatSetLocalToGlobalMapping() and MatSetStencil() first. The columns and rows in the stencil passed in MUST be contained within the ghost region of the given process as set with DACreateXXX() or MatSetStencil(). For example, if you create a DA with an overlap of one grid level and on a particular process its first local nonghost x logical coordinate is 6 (so its first ghost x logical coordinate is 5) the first i index you can use in your column and row indices in MatSetStencil() is 5. In Fortran idxm and idxn should be declared as $ MatStencil idxm(4,m),idxn(4,n) and the values inserted using $ idxm(MatStencil_i,1) = i $ idxm(MatStencil_j,1) = j $ idxm(MatStencil_k,1) = k etc Negative indices may be passed in idxm and idxn, these rows and columns are simply ignored. This allows easily inserting element stiffness matrices with homogeneous Dirchlet boundary conditions that you don't want represented in the matrix. Inspired by the structured grid interface to the HYPRE package (http://www.llnl.gov/CASC/hypre) Level: beginner Concepts: matrices^putting entries in .seealso: MatSetOption(), MatAssemblyBegin(), MatAssemblyEnd(), MatSetValuesBlocked(), MatSetValuesLocal() MatSetValues(), MatSetValuesStencil(), MatSetStencil(), DAGetMatrix(), DAVecGetArray(), MatStencil @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlockedStencil(Mat mat,PetscInt m,const MatStencil idxm[],PetscInt n,const MatStencil idxn[],const PetscScalar v[],InsertMode addv) { PetscErrorCode ierr; PetscInt j,i,jdxm[128],jdxn[256],dim = mat->stencil.dim,*dims = mat->stencil.dims+1,tmp; PetscInt *starts = mat->stencil.starts,*dxm = (PetscInt*)idxm,*dxn = (PetscInt*)idxn,sdim = dim - (1 - (PetscInt)mat->stencil.noc); PetscFunctionBegin; if (!m || !n) PetscFunctionReturn(0); /* no values to insert */ PetscValidHeaderSpecific(mat,MAT_COOKIE,1); PetscValidType(mat,1); PetscValidIntPointer(idxm,3); PetscValidIntPointer(idxn,5); PetscValidScalarPointer(v,6); if (m > 128) SETERRQ1(PETSC_ERR_SUP,"Can only set 128 rows at a time; trying to set %D",m); if (n > 128) SETERRQ1(PETSC_ERR_SUP,"Can only set 256 columns at a time; trying to set %D",n); for (i=0; istencil.dim = dim + (dof > 1); for (i=0; istencil.dims[i] = dims[dim-i-1]; /* copy the values in backwards */ mat->stencil.starts[i] = starts[dim-i-1]; } mat->stencil.dims[dim] = dof; mat->stencil.starts[dim] = 0; mat->stencil.noc = (PetscTruth)(dof == 1); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatSetValuesBlocked" /*@ MatSetValuesBlocked - Inserts or adds a block of values into a matrix. Not Collective Input Parameters: + mat - the matrix . v - a logically two-dimensional array of values . m, idxm - the number of block rows and their global block indices . n, idxn - the number of block columns and their global block indices - addv - either ADD_VALUES or INSERT_VALUES, where ADD_VALUES adds values to any existing entries, and INSERT_VALUES replaces existing entries with new values Notes: The m and n count the NUMBER of blocks in the row direction and column direction, NOT the total number of rows/columns; for example, if the block size is 2 and you are passing in values for rows 2,3,4,5 then m would be 2 (not 4). By default the values, v, are row-oriented and unsorted. So the layout of v is the same as for MatSetValues(). See MatSetOption() for other options. Calls to MatSetValuesBlocked() with the INSERT_VALUES and ADD_VALUES options cannot be mixed without intervening calls to the assembly routines. MatSetValuesBlocked() uses 0-based row and column numbers in Fortran as well as in C. Negative indices may be passed in idxm and idxn, these rows and columns are simply ignored. This allows easily inserting element stiffness matrices with homogeneous Dirchlet boundary conditions that you don't want represented in the matrix. Each time an entry is set within a sparse matrix via MatSetValues(), internal searching must be done to determine where to place the the data in the matrix storage space. By instead inserting blocks of entries via MatSetValuesBlocked(), the overhead of matrix assembly is reduced. Example: $ Suppose m=n=2 and block size(bs) = 2 The matrix is $ $ 1 2 | 3 4 $ 5 6 | 7 8 $ - - - | - - - $ 9 10 | 11 12 $ 13 14 | 15 16 $ $ v[] should be passed in like $ v[] = [1,2,3,4,5,6,7,8,9,10,11,12,13,14,15,16] Restrictions: MatSetValuesBlocked() is currently supported only for the BAIJ and SBAIJ formats Level: intermediate Concepts: matrices^putting entries in blocked .seealso: MatSetOption(), MatAssemblyBegin(), MatAssemblyEnd(), MatSetValues(), MatSetValuesBlockedLocal() @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlocked(Mat mat,PetscInt m,const PetscInt idxm[],PetscInt n,const PetscInt idxn[],const PetscScalar v[],InsertMode addv) { PetscErrorCode ierr; PetscFunctionBegin; if (!m || !n) PetscFunctionReturn(0); /* no values to insert */ PetscValidHeaderSpecific(mat,MAT_COOKIE,1); PetscValidType(mat,1); PetscValidIntPointer(idxm,3); PetscValidIntPointer(idxn,5); PetscValidScalarPointer(v,6); ierr = MatPreallocated(mat);CHKERRQ(ierr); if (mat->insertmode == NOT_SET_VALUES) { mat->insertmode = addv; } #if defined(PETSC_USE_DEBUG) else if (mat->insertmode != addv) { SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Cannot mix add values and insert values"); } if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); #endif if (mat->assembled) { mat->was_assembled = PETSC_TRUE; mat->assembled = PETSC_FALSE; } ierr = PetscLogEventBegin(MAT_SetValues,mat,0,0,0);CHKERRQ(ierr); if (!mat->ops->setvaluesblocked) SETERRQ1(PETSC_ERR_SUP,"Mat type %s",mat->type_name); ierr = (*mat->ops->setvaluesblocked)(mat,m,idxm,n,idxn,v,addv);CHKERRQ(ierr); ierr = PetscLogEventEnd(MAT_SetValues,mat,0,0,0);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatGetValues" /*@ MatGetValues - Gets a block of values from a matrix. Not Collective; currently only returns a local block Input Parameters: + mat - the matrix . v - a logically two-dimensional array for storing the values . m, idxm - the number of rows and their global indices - n, idxn - the number of columns and their global indices Notes: The user must allocate space (m*n PetscScalars) for the values, v. The values, v, are then returned in a row-oriented format, analogous to that used by default in MatSetValues(). MatGetValues() uses 0-based row and column numbers in Fortran as well as in C. MatGetValues() requires that the matrix has been assembled with MatAssemblyBegin()/MatAssemblyEnd(). Thus, calls to MatSetValues() and MatGetValues() CANNOT be made in succession without intermediate matrix assembly. Level: advanced Concepts: matrices^accessing values .seealso: MatGetRow(), MatGetSubMatrices(), MatSetValues() @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatGetValues(Mat mat,PetscInt m,const PetscInt idxm[],PetscInt n,const PetscInt idxn[],PetscScalar v[]) { PetscErrorCode ierr; PetscFunctionBegin; PetscValidHeaderSpecific(mat,MAT_COOKIE,1); PetscValidType(mat,1); PetscValidIntPointer(idxm,3); PetscValidIntPointer(idxn,5); PetscValidScalarPointer(v,6); if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); if (!mat->ops->getvalues) SETERRQ1(PETSC_ERR_SUP,"Mat type %s",mat->type_name); ierr = MatPreallocated(mat);CHKERRQ(ierr); ierr = PetscLogEventBegin(MAT_GetValues,mat,0,0,0);CHKERRQ(ierr); ierr = (*mat->ops->getvalues)(mat,m,idxm,n,idxn,v);CHKERRQ(ierr); ierr = PetscLogEventEnd(MAT_GetValues,mat,0,0,0);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatSetLocalToGlobalMapping" /*@ MatSetLocalToGlobalMapping - Sets a local-to-global numbering for use by the routine MatSetValuesLocal() to allow users to insert matrix entries using a local (per-processor) numbering. Not Collective Input Parameters: + x - the matrix - mapping - mapping created with ISLocalToGlobalMappingCreate() or ISLocalToGlobalMappingCreateIS() Level: intermediate Concepts: matrices^local to global mapping Concepts: local to global mapping^for matrices .seealso: MatAssemblyBegin(), MatAssemblyEnd(), MatSetValues(), MatSetValuesLocal() @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatSetLocalToGlobalMapping(Mat x,ISLocalToGlobalMapping mapping) { PetscErrorCode ierr; PetscFunctionBegin; PetscValidHeaderSpecific(x,MAT_COOKIE,1); PetscValidType(x,1); PetscValidHeaderSpecific(mapping,IS_LTOGM_COOKIE,2); if (x->mapping) { SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Mapping already set for matrix"); } ierr = MatPreallocated(x);CHKERRQ(ierr); if (x->ops->setlocaltoglobalmapping) { ierr = (*x->ops->setlocaltoglobalmapping)(x,mapping);CHKERRQ(ierr); } else { ierr = PetscObjectReference((PetscObject)mapping);CHKERRQ(ierr); if (x->mapping) { ierr = ISLocalToGlobalMappingDestroy(x->mapping);CHKERRQ(ierr); } x->mapping = mapping; } PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatSetLocalToGlobalMappingBlock" /*@ MatSetLocalToGlobalMappingBlock - Sets a local-to-global numbering for use by the routine MatSetValuesBlockedLocal() to allow users to insert matrix entries using a local (per-processor) numbering. Not Collective Input Parameters: + x - the matrix - mapping - mapping created with ISLocalToGlobalMappingCreate() or ISLocalToGlobalMappingCreateIS() Level: intermediate Concepts: matrices^local to global mapping blocked Concepts: local to global mapping^for matrices, blocked .seealso: MatAssemblyBegin(), MatAssemblyEnd(), MatSetValues(), MatSetValuesBlockedLocal(), MatSetValuesBlocked(), MatSetValuesLocal() @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatSetLocalToGlobalMappingBlock(Mat x,ISLocalToGlobalMapping mapping) { PetscErrorCode ierr; PetscFunctionBegin; PetscValidHeaderSpecific(x,MAT_COOKIE,1); PetscValidType(x,1); PetscValidHeaderSpecific(mapping,IS_LTOGM_COOKIE,2); if (x->bmapping) { SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Mapping already set for matrix"); } ierr = PetscObjectReference((PetscObject)mapping);CHKERRQ(ierr); if (x->bmapping) { ierr = ISLocalToGlobalMappingDestroy(x->mapping);CHKERRQ(ierr); } x->bmapping = mapping; PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatSetValuesLocal" /*@ MatSetValuesLocal - Inserts or adds values into certain locations of a matrix, using a local ordering of the nodes. Not Collective Input Parameters: + x - the matrix . nrow, irow - number of rows and their local indices . ncol, icol - number of columns and their local indices . y - a logically two-dimensional array of values - addv - either INSERT_VALUES or ADD_VALUES, where ADD_VALUES adds values to any existing entries, and INSERT_VALUES replaces existing entries with new values Notes: Before calling MatSetValuesLocal(), the user must first set the local-to-global mapping by calling MatSetLocalToGlobalMapping(). Calls to MatSetValuesLocal() with the INSERT_VALUES and ADD_VALUES options cannot be mixed without intervening calls to the assembly routines. These values may be cached, so MatAssemblyBegin() and MatAssemblyEnd() MUST be called after all calls to MatSetValuesLocal() have been completed. Level: intermediate Concepts: matrices^putting entries in with local numbering .seealso: MatAssemblyBegin(), MatAssemblyEnd(), MatSetValues(), MatSetLocalToGlobalMapping(), MatSetValueLocal() @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesLocal(Mat mat,PetscInt nrow,const PetscInt irow[],PetscInt ncol,const PetscInt icol[],const PetscScalar y[],InsertMode addv) { PetscErrorCode ierr; PetscInt irowm[2048],icolm[2048]; PetscFunctionBegin; PetscValidHeaderSpecific(mat,MAT_COOKIE,1); PetscValidType(mat,1); PetscValidIntPointer(irow,3); PetscValidIntPointer(icol,5); PetscValidScalarPointer(y,6); ierr = MatPreallocated(mat);CHKERRQ(ierr); if (mat->insertmode == NOT_SET_VALUES) { mat->insertmode = addv; } #if defined(PETSC_USE_DEBUG) else if (mat->insertmode != addv) { SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Cannot mix add values and insert values"); } if (!mat->ops->setvalueslocal && (nrow > 2048 || ncol > 2048)) { SETERRQ2(PETSC_ERR_SUP,"Number column/row indices must be <= 2048: are %D %D",nrow,ncol); } if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); #endif if (mat->assembled) { mat->was_assembled = PETSC_TRUE; mat->assembled = PETSC_FALSE; } ierr = PetscLogEventBegin(MAT_SetValues,mat,0,0,0);CHKERRQ(ierr); if (!mat->ops->setvalueslocal) { ierr = ISLocalToGlobalMappingApply(mat->mapping,nrow,irow,irowm);CHKERRQ(ierr); ierr = ISLocalToGlobalMappingApply(mat->mapping,ncol,icol,icolm);CHKERRQ(ierr); ierr = (*mat->ops->setvalues)(mat,nrow,irowm,ncol,icolm,y,addv);CHKERRQ(ierr); } else { ierr = (*mat->ops->setvalueslocal)(mat,nrow,irow,ncol,icol,y,addv);CHKERRQ(ierr); } mat->same_nonzero = PETSC_FALSE; ierr = PetscLogEventEnd(MAT_SetValues,mat,0,0,0);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatSetValuesBlockedLocal" /*@ MatSetValuesBlockedLocal - Inserts or adds values into certain locations of a matrix, using a local ordering of the nodes a block at a time. Not Collective Input Parameters: + x - the matrix . nrow, irow - number of rows and their local indices . ncol, icol - number of columns and their local indices . y - a logically two-dimensional array of values - addv - either INSERT_VALUES or ADD_VALUES, where ADD_VALUES adds values to any existing entries, and INSERT_VALUES replaces existing entries with new values Notes: Before calling MatSetValuesBlockedLocal(), the user must first set the local-to-global mapping by calling MatSetLocalToGlobalMappingBlock(), where the mapping MUST be set for matrix blocks, not for matrix elements. Calls to MatSetValuesBlockedLocal() with the INSERT_VALUES and ADD_VALUES options cannot be mixed without intervening calls to the assembly routines. These values may be cached, so MatAssemblyBegin() and MatAssemblyEnd() MUST be called after all calls to MatSetValuesBlockedLocal() have been completed. Level: intermediate Concepts: matrices^putting blocked values in with local numbering .seealso: MatAssemblyBegin(), MatAssemblyEnd(), MatSetValuesLocal(), MatSetLocalToGlobalMappingBlock(), MatSetValuesBlocked() @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlockedLocal(Mat mat,PetscInt nrow,const PetscInt irow[],PetscInt ncol,const PetscInt icol[],const PetscScalar y[],InsertMode addv) { PetscErrorCode ierr; PetscInt irowm[2048],icolm[2048]; PetscFunctionBegin; PetscValidHeaderSpecific(mat,MAT_COOKIE,1); PetscValidType(mat,1); PetscValidIntPointer(irow,3); PetscValidIntPointer(icol,5); PetscValidScalarPointer(y,6); ierr = MatPreallocated(mat);CHKERRQ(ierr); if (mat->insertmode == NOT_SET_VALUES) { mat->insertmode = addv; } #if defined(PETSC_USE_DEBUG) else if (mat->insertmode != addv) { SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Cannot mix add values and insert values"); } if (!mat->bmapping) { SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Local to global never set with MatSetLocalToGlobalMappingBlock()"); } if (nrow > 2048 || ncol > 2048) { SETERRQ2(PETSC_ERR_SUP,"Number column/row indices must be <= 2048: are %D %D",nrow,ncol); } if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); #endif if (mat->assembled) { mat->was_assembled = PETSC_TRUE; mat->assembled = PETSC_FALSE; } ierr = PetscLogEventBegin(MAT_SetValues,mat,0,0,0);CHKERRQ(ierr); ierr = ISLocalToGlobalMappingApply(mat->bmapping,nrow,irow,irowm);CHKERRQ(ierr); ierr = ISLocalToGlobalMappingApply(mat->bmapping,ncol,icol,icolm);CHKERRQ(ierr); ierr = (*mat->ops->setvaluesblocked)(mat,nrow,irowm,ncol,icolm,y,addv);CHKERRQ(ierr); ierr = PetscLogEventEnd(MAT_SetValues,mat,0,0,0);CHKERRQ(ierr); PetscFunctionReturn(0); } /* --------------------------------------------------------*/ #undef __FUNCT__ #define __FUNCT__ "MatMult" /*@ MatMult - Computes the matrix-vector product, y = Ax. Collective on Mat and Vec Input Parameters: + mat - the matrix - x - the vector to be multiplied Output Parameters: . y - the result Notes: The vectors x and y cannot be the same. I.e., one cannot call MatMult(A,y,y). Level: beginner Concepts: matrix-vector product .seealso: MatMultTranspose(), MatMultAdd(), MatMultTransposeAdd() @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatMult(Mat mat,Vec x,Vec y) { PetscErrorCode ierr; PetscFunctionBegin; PetscValidHeaderSpecific(mat,MAT_COOKIE,1); PetscValidType(mat,1); PetscValidHeaderSpecific(x,VEC_COOKIE,2); PetscValidHeaderSpecific(y,VEC_COOKIE,3); if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); if (x == y) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"x and y must be different vectors"); #ifndef PETSC_HAVE_CONSTRAINTS if (mat->cmap.N != x->map.N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec x: global dim %D %D",mat->cmap.N,x->map.N); if (mat->rmap.N != y->map.N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec y: global dim %D %D",mat->rmap.N,y->map.N); if (mat->rmap.n != y->map.n) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec y: local dim %D %D",mat->rmap.n,y->map.n); #endif ierr = MatPreallocated(mat);CHKERRQ(ierr); if (mat->nullsp) { ierr = MatNullSpaceRemove(mat->nullsp,x,&x);CHKERRQ(ierr); } if (!mat->ops->mult) SETERRQ(PETSC_ERR_SUP,"This matrix type does not have a multiply defined"); ierr = PetscLogEventBegin(MAT_Mult,mat,x,y,0);CHKERRQ(ierr); ierr = (*mat->ops->mult)(mat,x,y);CHKERRQ(ierr); ierr = PetscLogEventEnd(MAT_Mult,mat,x,y,0);CHKERRQ(ierr); if (mat->nullsp) { ierr = MatNullSpaceRemove(mat->nullsp,y,PETSC_NULL);CHKERRQ(ierr); } ierr = PetscObjectStateIncrease((PetscObject)y);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatMultTranspose" /*@ MatMultTranspose - Computes matrix transpose times a vector. Collective on Mat and Vec Input Parameters: + mat - the matrix - x - the vector to be multilplied Output Parameters: . y - the result Notes: The vectors x and y cannot be the same. I.e., one cannot call MatMultTranspose(A,y,y). Level: beginner Concepts: matrix vector product^transpose .seealso: MatMult(), MatMultAdd(), MatMultTransposeAdd() @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatMultTranspose(Mat mat,Vec x,Vec y) { PetscErrorCode ierr; PetscFunctionBegin; PetscValidHeaderSpecific(mat,MAT_COOKIE,1); PetscValidType(mat,1); PetscValidHeaderSpecific(x,VEC_COOKIE,2); PetscValidHeaderSpecific(y,VEC_COOKIE,3); if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); if (x == y) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"x and y must be different vectors"); #ifndef PETSC_HAVE_CONSTRAINTS if (mat->rmap.N != x->map.N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec x: global dim %D %D",mat->rmap.N,x->map.N); if (mat->cmap.N != y->map.N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec y: global dim %D %D",mat->cmap.N,y->map.N); #endif ierr = MatPreallocated(mat);CHKERRQ(ierr); if (!mat->ops->multtranspose) SETERRQ(PETSC_ERR_SUP,"This matrix type does not have a multiply tranpose defined"); ierr = PetscLogEventBegin(MAT_MultTranspose,mat,x,y,0);CHKERRQ(ierr); ierr = (*mat->ops->multtranspose)(mat,x,y);CHKERRQ(ierr); ierr = PetscLogEventEnd(MAT_MultTranspose,mat,x,y,0);CHKERRQ(ierr); ierr = PetscObjectStateIncrease((PetscObject)y);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatMultAdd" /*@ MatMultAdd - Computes v3 = v2 + A * v1. Collective on Mat and Vec Input Parameters: + mat - the matrix - v1, v2 - the vectors Output Parameters: . v3 - the result Notes: The vectors v1 and v3 cannot be the same. I.e., one cannot call MatMultAdd(A,v1,v2,v1). Level: beginner Concepts: matrix vector product^addition .seealso: MatMultTranspose(), MatMult(), MatMultTransposeAdd() @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatMultAdd(Mat mat,Vec v1,Vec v2,Vec v3) { PetscErrorCode ierr; PetscFunctionBegin; PetscValidHeaderSpecific(mat,MAT_COOKIE,1); PetscValidType(mat,1); PetscValidHeaderSpecific(v1,VEC_COOKIE,2); PetscValidHeaderSpecific(v2,VEC_COOKIE,3); PetscValidHeaderSpecific(v3,VEC_COOKIE,4); if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); if (mat->cmap.N != v1->map.N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec v1: global dim %D %D",mat->cmap.N,v1->map.N); if (mat->rmap.N != v2->map.N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec v2: global dim %D %D",mat->rmap.N,v2->map.N); if (mat->rmap.N != v3->map.N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec v3: global dim %D %D",mat->rmap.N,v3->map.N); if (mat->rmap.n != v3->map.n) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec v3: local dim %D %D",mat->rmap.n,v3->map.n); if (mat->rmap.n != v2->map.n) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec v2: local dim %D %D",mat->rmap.n,v2->map.n); if (v1 == v3) SETERRQ(PETSC_ERR_ARG_IDN,"v1 and v3 must be different vectors"); ierr = MatPreallocated(mat);CHKERRQ(ierr); ierr = PetscLogEventBegin(MAT_MultAdd,mat,v1,v2,v3);CHKERRQ(ierr); ierr = (*mat->ops->multadd)(mat,v1,v2,v3);CHKERRQ(ierr); ierr = PetscLogEventEnd(MAT_MultAdd,mat,v1,v2,v3);CHKERRQ(ierr); ierr = PetscObjectStateIncrease((PetscObject)v3);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatMultTransposeAdd" /*@ MatMultTransposeAdd - Computes v3 = v2 + A' * v1. Collective on Mat and Vec Input Parameters: + mat - the matrix - v1, v2 - the vectors Output Parameters: . v3 - the result Notes: The vectors v1 and v3 cannot be the same. I.e., one cannot call MatMultTransposeAdd(A,v1,v2,v1). Level: beginner Concepts: matrix vector product^transpose and addition .seealso: MatMultTranspose(), MatMultAdd(), MatMult() @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeAdd(Mat mat,Vec v1,Vec v2,Vec v3) { PetscErrorCode ierr; PetscFunctionBegin; PetscValidHeaderSpecific(mat,MAT_COOKIE,1); PetscValidType(mat,1); PetscValidHeaderSpecific(v1,VEC_COOKIE,2); PetscValidHeaderSpecific(v2,VEC_COOKIE,3); PetscValidHeaderSpecific(v3,VEC_COOKIE,4); if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); if (!mat->ops->multtransposeadd) SETERRQ1(PETSC_ERR_SUP,"Mat type %s",mat->type_name); if (v1 == v3) SETERRQ(PETSC_ERR_ARG_IDN,"v1 and v3 must be different vectors"); if (mat->rmap.N != v1->map.N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec v1: global dim %D %D",mat->rmap.N,v1->map.N); if (mat->cmap.N != v2->map.N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec v2: global dim %D %D",mat->cmap.N,v2->map.N); if (mat->cmap.N != v3->map.N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec v3: global dim %D %D",mat->cmap.N,v3->map.N); ierr = MatPreallocated(mat);CHKERRQ(ierr); ierr = PetscLogEventBegin(MAT_MultTransposeAdd,mat,v1,v2,v3);CHKERRQ(ierr); ierr = (*mat->ops->multtransposeadd)(mat,v1,v2,v3);CHKERRQ(ierr); ierr = PetscLogEventEnd(MAT_MultTransposeAdd,mat,v1,v2,v3);CHKERRQ(ierr); ierr = PetscObjectStateIncrease((PetscObject)v3);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatMultConstrained" /*@ MatMultConstrained - The inner multiplication routine for a constrained matrix P^T A P. Collective on Mat and Vec Input Parameters: + mat - the matrix - x - the vector to be multilplied Output Parameters: . y - the result Notes: The vectors x and y cannot be the same. I.e., one cannot call MatMult(A,y,y). Level: beginner .keywords: matrix, multiply, matrix-vector product, constraint .seealso: MatMult(), MatMultTrans(), MatMultAdd(), MatMultTransAdd() @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatMultConstrained(Mat mat,Vec x,Vec y) { PetscErrorCode ierr; PetscFunctionBegin; PetscValidHeaderSpecific(mat,MAT_COOKIE,1); PetscValidHeaderSpecific(x,VEC_COOKIE,2); PetscValidHeaderSpecific(y,VEC_COOKIE,3); if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); if (x == y) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"x and y must be different vectors"); if (mat->cmap.N != x->map.N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec x: global dim %D %D",mat->cmap.N,x->map.N); if (mat->rmap.N != y->map.N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec y: global dim %D %D",mat->rmap.N,y->map.N); if (mat->rmap.n != y->map.n) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec y: local dim %D %D",mat->rmap.n,y->map.n); ierr = PetscLogEventBegin(MAT_MultConstrained,mat,x,y,0);CHKERRQ(ierr); ierr = (*mat->ops->multconstrained)(mat,x,y);CHKERRQ(ierr); ierr = PetscLogEventEnd(MAT_MultConstrained,mat,x,y,0);CHKERRQ(ierr); ierr = PetscObjectStateIncrease((PetscObject)y);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatMultTransposeConstrained" /*@ MatMultTransposeConstrained - The inner multiplication routine for a constrained matrix P^T A^T P. Collective on Mat and Vec Input Parameters: + mat - the matrix - x - the vector to be multilplied Output Parameters: . y - the result Notes: The vectors x and y cannot be the same. I.e., one cannot call MatMult(A,y,y). Level: beginner .keywords: matrix, multiply, matrix-vector product, constraint .seealso: MatMult(), MatMultTrans(), MatMultAdd(), MatMultTransAdd() @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeConstrained(Mat mat,Vec x,Vec y) { PetscErrorCode ierr; PetscFunctionBegin; PetscValidHeaderSpecific(mat,MAT_COOKIE,1); PetscValidHeaderSpecific(x,VEC_COOKIE,2); PetscValidHeaderSpecific(y,VEC_COOKIE,3); if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); if (x == y) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"x and y must be different vectors"); if (mat->rmap.N != x->map.N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec x: global dim %D %D",mat->cmap.N,x->map.N); if (mat->cmap.N != y->map.N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec y: global dim %D %D",mat->rmap.N,y->map.N); ierr = PetscLogEventBegin(MAT_MultConstrained,mat,x,y,0);CHKERRQ(ierr); ierr = (*mat->ops->multtransposeconstrained)(mat,x,y);CHKERRQ(ierr); ierr = PetscLogEventEnd(MAT_MultConstrained,mat,x,y,0);CHKERRQ(ierr); ierr = PetscObjectStateIncrease((PetscObject)y);CHKERRQ(ierr); PetscFunctionReturn(0); } /* ------------------------------------------------------------*/ #undef __FUNCT__ #define __FUNCT__ "MatGetInfo" /*@ MatGetInfo - Returns information about matrix storage (number of nonzeros, memory, etc.). Collective on Mat if MAT_GLOBAL_MAX or MAT_GLOBAL_SUM is used as the flag Input Parameters: . mat - the matrix Output Parameters: + flag - flag indicating the type of parameters to be returned (MAT_LOCAL - local matrix, MAT_GLOBAL_MAX - maximum over all processors, MAT_GLOBAL_SUM - sum over all processors) - info - matrix information context Notes: The MatInfo context contains a variety of matrix data, including number of nonzeros allocated and used, number of mallocs during matrix assembly, etc. Additional information for factored matrices is provided (such as the fill ratio, number of mallocs during factorization, etc.). Much of this info is printed to STDOUT when using the runtime options $ -info -mat_view_info Example for C/C++ Users: See the file ${PETSC_DIR}/include/petscmat.h for a complete list of data within the MatInfo context. For example, .vb MatInfo info; Mat A; double mal, nz_a, nz_u; MatGetInfo(A,MAT_LOCAL,&info); mal = info.mallocs; nz_a = info.nz_allocated; .ve Example for Fortran Users: Fortran users should declare info as a double precision array of dimension MAT_INFO_SIZE, and then extract the parameters of interest. See the file ${PETSC_DIR}/include/finclude/petscmat.h a complete list of parameter names. .vb double precision info(MAT_INFO_SIZE) double precision mal, nz_a Mat A integer ierr call MatGetInfo(A,MAT_LOCAL,info,ierr) mal = info(MAT_INFO_MALLOCS) nz_a = info(MAT_INFO_NZ_ALLOCATED) .ve Level: intermediate Concepts: matrices^getting information on @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatGetInfo(Mat mat,MatInfoType flag,MatInfo *info) { PetscErrorCode ierr; PetscFunctionBegin; PetscValidHeaderSpecific(mat,MAT_COOKIE,1); PetscValidType(mat,1); PetscValidPointer(info,3); if (!mat->ops->getinfo) SETERRQ1(PETSC_ERR_SUP,"Mat type %s",mat->type_name); ierr = MatPreallocated(mat);CHKERRQ(ierr); ierr = (*mat->ops->getinfo)(mat,flag,info);CHKERRQ(ierr); PetscFunctionReturn(0); } /* ----------------------------------------------------------*/ #undef __FUNCT__ #define __FUNCT__ "MatILUDTFactor" /*@C MatILUDTFactor - Performs a drop tolerance ILU factorization. Collective on Mat Input Parameters: + mat - the matrix . row - row permutation . col - column permutation - info - information about the factorization to be done Output Parameters: . fact - the factored matrix Level: developer Notes: Most users should employ the simplified KSP interface for linear solvers instead of working directly with matrix algebra routines such as this. See, e.g., KSPCreate(). This is currently only supported for the SeqAIJ matrix format using code from Yousef Saad's SPARSEKIT2 package (translated to C with f2c) and/or Matlab. SPARSEKIT2 is copyrighted by Yousef Saad with the GNU copyright and thus can be distributed with PETSc. Concepts: matrices^ILUDT factorization .seealso: MatLUFactorSymbolic(), MatLUFactorNumeric(), MatCholeskyFactor(), MatFactorInfo @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatILUDTFactor(Mat mat,IS row,IS col,MatFactorInfo *info,Mat *fact) { PetscErrorCode ierr; PetscFunctionBegin; PetscValidHeaderSpecific(mat,MAT_COOKIE,1); PetscValidType(mat,1); if (row) PetscValidHeaderSpecific(row,IS_COOKIE,2); if (col) PetscValidHeaderSpecific(col,IS_COOKIE,3); PetscValidPointer(info,4); PetscValidPointer(fact,5); if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); if (!mat->ops->iludtfactor) SETERRQ1(PETSC_ERR_SUP,"Mat type %s",mat->type_name); ierr = MatPreallocated(mat);CHKERRQ(ierr); ierr = PetscLogEventBegin(MAT_ILUFactor,mat,row,col,0);CHKERRQ(ierr); ierr = (*mat->ops->iludtfactor)(mat,row,col,info,fact);CHKERRQ(ierr); ierr = PetscLogEventEnd(MAT_ILUFactor,mat,row,col,0);CHKERRQ(ierr); ierr = PetscObjectStateIncrease((PetscObject)*fact);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatLUFactor" /*@ MatLUFactor - Performs in-place LU factorization of matrix. Collective on Mat Input Parameters: + mat - the matrix . row - row permutation . col - column permutation - info - options for factorization, includes $ fill - expected fill as ratio of original fill. $ dtcol - pivot tolerance (0 no pivot, 1 full column pivoting) $ Run with the option -info to determine an optimal value to use Notes: Most users should employ the simplified KSP interface for linear solvers instead of working directly with matrix algebra routines such as this. See, e.g., KSPCreate(). This changes the state of the matrix to a factored matrix; it cannot be used for example with MatSetValues() unless one first calls MatSetUnfactored(). Level: developer Concepts: matrices^LU factorization .seealso: MatLUFactorSymbolic(), MatLUFactorNumeric(), MatCholeskyFactor(), MatGetOrdering(), MatSetUnfactored(), MatFactorInfo @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactor(Mat mat,IS row,IS col,MatFactorInfo *info) { PetscErrorCode ierr; PetscFunctionBegin; PetscValidHeaderSpecific(mat,MAT_COOKIE,1); if (row) PetscValidHeaderSpecific(row,IS_COOKIE,2); if (col) PetscValidHeaderSpecific(col,IS_COOKIE,3); PetscValidPointer(info,4); PetscValidType(mat,1); if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); if (!mat->ops->lufactor) SETERRQ1(PETSC_ERR_SUP,"Mat type %s",mat->type_name); ierr = MatPreallocated(mat);CHKERRQ(ierr); ierr = PetscLogEventBegin(MAT_LUFactor,mat,row,col,0);CHKERRQ(ierr); ierr = (*mat->ops->lufactor)(mat,row,col,info);CHKERRQ(ierr); ierr = PetscLogEventEnd(MAT_LUFactor,mat,row,col,0);CHKERRQ(ierr); ierr = PetscObjectStateIncrease((PetscObject)mat);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatILUFactor" /*@ MatILUFactor - Performs in-place ILU factorization of matrix. Collective on Mat Input Parameters: + mat - the matrix . row - row permutation . col - column permutation - info - structure containing $ levels - number of levels of fill. $ expected fill - as ratio of original fill. $ 1 or 0 - indicating force fill on diagonal (improves robustness for matrices missing diagonal entries) Notes: Probably really in-place only when level of fill is zero, otherwise allocates new space to store factored matrix and deletes previous memory. Most users should employ the simplified KSP interface for linear solvers instead of working directly with matrix algebra routines such as this. See, e.g., KSPCreate(). Level: developer Concepts: matrices^ILU factorization .seealso: MatILUFactorSymbolic(), MatLUFactorNumeric(), MatCholeskyFactor(), MatFactorInfo @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatILUFactor(Mat mat,IS row,IS col,MatFactorInfo *info) { PetscErrorCode ierr; PetscFunctionBegin; PetscValidHeaderSpecific(mat,MAT_COOKIE,1); if (row) PetscValidHeaderSpecific(row,IS_COOKIE,2); if (col) PetscValidHeaderSpecific(col,IS_COOKIE,3); PetscValidPointer(info,4); PetscValidType(mat,1); if (mat->rmap.N != mat->cmap.N) SETERRQ(PETSC_ERR_ARG_WRONG,"matrix must be square"); if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); if (!mat->ops->ilufactor) SETERRQ1(PETSC_ERR_SUP,"Mat type %s",mat->type_name); ierr = MatPreallocated(mat);CHKERRQ(ierr); ierr = PetscLogEventBegin(MAT_ILUFactor,mat,row,col,0);CHKERRQ(ierr); ierr = (*mat->ops->ilufactor)(mat,row,col,info);CHKERRQ(ierr); ierr = PetscLogEventEnd(MAT_ILUFactor,mat,row,col,0);CHKERRQ(ierr); ierr = PetscObjectStateIncrease((PetscObject)mat);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatLUFactorSymbolic" /*@ MatLUFactorSymbolic - Performs symbolic LU factorization of matrix. Call this routine before calling MatLUFactorNumeric(). Collective on Mat Input Parameters: + mat - the matrix . row, col - row and column permutations - info - options for factorization, includes $ fill - expected fill as ratio of original fill. $ dtcol - pivot tolerance (0 no pivot, 1 full column pivoting) $ Run with the option -info to determine an optimal value to use Output Parameter: . fact - new matrix that has been symbolically factored Notes: See the users manual for additional information about choosing the fill factor for better efficiency. Most users should employ the simplified KSP interface for linear solvers instead of working directly with matrix algebra routines such as this. See, e.g., KSPCreate(). Level: developer Concepts: matrices^LU symbolic factorization .seealso: MatLUFactor(), MatLUFactorNumeric(), MatCholeskyFactor(), MatFactorInfo @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactorSymbolic(Mat mat,IS row,IS col,MatFactorInfo *info,Mat *fact) { PetscErrorCode ierr; PetscFunctionBegin; PetscValidHeaderSpecific(mat,MAT_COOKIE,1); if (row) PetscValidHeaderSpecific(row,IS_COOKIE,2); if (col) PetscValidHeaderSpecific(col,IS_COOKIE,3); PetscValidPointer(info,4); PetscValidType(mat,1); PetscValidPointer(fact,5); if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); if (!mat->ops->lufactorsymbolic) SETERRQ1(PETSC_ERR_SUP,"Matrix type %s symbolic LU",mat->type_name); ierr = MatPreallocated(mat);CHKERRQ(ierr); ierr = PetscLogEventBegin(MAT_LUFactorSymbolic,mat,row,col,0);CHKERRQ(ierr); ierr = (*mat->ops->lufactorsymbolic)(mat,row,col,info,fact);CHKERRQ(ierr); ierr = PetscLogEventEnd(MAT_LUFactorSymbolic,mat,row,col,0);CHKERRQ(ierr); ierr = PetscObjectStateIncrease((PetscObject)*fact);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatLUFactorNumeric" /*@ MatLUFactorNumeric - Performs numeric LU factorization of a matrix. Call this routine after first calling MatLUFactorSymbolic(). Collective on Mat Input Parameters: + mat - the matrix . info - options for factorization - fact - the matrix generated for the factor, from MatLUFactorSymbolic() Notes: See MatLUFactor() for in-place factorization. See MatCholeskyFactorNumeric() for the symmetric, positive definite case. Most users should employ the simplified KSP interface for linear solvers instead of working directly with matrix algebra routines such as this. See, e.g., KSPCreate(). Level: developer Concepts: matrices^LU numeric factorization .seealso: MatLUFactorSymbolic(), MatLUFactor(), MatCholeskyFactor() @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactorNumeric(Mat mat,MatFactorInfo *info,Mat *fact) { PetscErrorCode ierr; PetscFunctionBegin; PetscValidHeaderSpecific(mat,MAT_COOKIE,1); PetscValidType(mat,1); PetscValidPointer(fact,2); PetscValidHeaderSpecific(*fact,MAT_COOKIE,2); if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); if (mat->rmap.N != (*fact)->rmap.N || mat->cmap.N != (*fact)->cmap.N) { SETERRQ4(PETSC_ERR_ARG_SIZ,"Mat mat,Mat *fact: global dimensions are different %D should = %D %D should = %D",mat->rmap.N,(*fact)->rmap.N,mat->cmap.N,(*fact)->cmap.N); } if (!(*fact)->ops->lufactornumeric) SETERRQ1(PETSC_ERR_SUP,"Mat type %s",mat->type_name); ierr = MatPreallocated(mat);CHKERRQ(ierr); ierr = PetscLogEventBegin(MAT_LUFactorNumeric,mat,*fact,0,0);CHKERRQ(ierr); ierr = (*(*fact)->ops->lufactornumeric)(mat,info,fact);CHKERRQ(ierr); ierr = PetscLogEventEnd(MAT_LUFactorNumeric,mat,*fact,0,0);CHKERRQ(ierr); ierr = MatView_Private(*fact);CHKERRQ(ierr); ierr = PetscObjectStateIncrease((PetscObject)*fact);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatCholeskyFactor" /*@ MatCholeskyFactor - Performs in-place Cholesky factorization of a symmetric matrix. Collective on Mat Input Parameters: + mat - the matrix . perm - row and column permutations - f - expected fill as ratio of original fill Notes: See MatLUFactor() for the nonsymmetric case. See also MatCholeskyFactorSymbolic(), and MatCholeskyFactorNumeric(). Most users should employ the simplified KSP interface for linear solvers instead of working directly with matrix algebra routines such as this. See, e.g., KSPCreate(). Level: developer Concepts: matrices^Cholesky factorization .seealso: MatLUFactor(), MatCholeskyFactorSymbolic(), MatCholeskyFactorNumeric() MatGetOrdering() @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactor(Mat mat,IS perm,MatFactorInfo *info) { PetscErrorCode ierr; PetscFunctionBegin; PetscValidHeaderSpecific(mat,MAT_COOKIE,1); PetscValidType(mat,1); PetscValidHeaderSpecific(perm,IS_COOKIE,2); PetscValidPointer(info,3); if (mat->rmap.N != mat->cmap.N) SETERRQ(PETSC_ERR_ARG_WRONG,"Matrix must be square"); if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); if (!mat->ops->choleskyfactor) SETERRQ1(PETSC_ERR_SUP,"Mat type %s",mat->type_name); ierr = MatPreallocated(mat);CHKERRQ(ierr); ierr = PetscLogEventBegin(MAT_CholeskyFactor,mat,perm,0,0);CHKERRQ(ierr); ierr = (*mat->ops->choleskyfactor)(mat,perm,info);CHKERRQ(ierr); ierr = PetscLogEventEnd(MAT_CholeskyFactor,mat,perm,0,0);CHKERRQ(ierr); ierr = PetscObjectStateIncrease((PetscObject)mat);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatCholeskyFactorSymbolic" /*@ MatCholeskyFactorSymbolic - Performs symbolic Cholesky factorization of a symmetric matrix. Collective on Mat Input Parameters: + mat - the matrix . perm - row and column permutations - info - options for factorization, includes $ fill - expected fill as ratio of original fill. $ dtcol - pivot tolerance (0 no pivot, 1 full column pivoting) $ Run with the option -info to determine an optimal value to use Output Parameter: . fact - the factored matrix Notes: See MatLUFactorSymbolic() for the nonsymmetric case. See also MatCholeskyFactor() and MatCholeskyFactorNumeric(). Most users should employ the simplified KSP interface for linear solvers instead of working directly with matrix algebra routines such as this. See, e.g., KSPCreate(). Level: developer Concepts: matrices^Cholesky symbolic factorization .seealso: MatLUFactorSymbolic(), MatCholeskyFactor(), MatCholeskyFactorNumeric() MatGetOrdering() @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactorSymbolic(Mat mat,IS perm,MatFactorInfo *info,Mat *fact) { PetscErrorCode ierr; PetscFunctionBegin; PetscValidHeaderSpecific(mat,MAT_COOKIE,1); PetscValidType(mat,1); if (perm) PetscValidHeaderSpecific(perm,IS_COOKIE,2); PetscValidPointer(info,3); PetscValidPointer(fact,4); if (mat->rmap.N != mat->cmap.N) SETERRQ(PETSC_ERR_ARG_WRONG,"Matrix must be square"); if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); if (!mat->ops->choleskyfactorsymbolic) SETERRQ1(PETSC_ERR_SUP,"Mat type %s",mat->type_name); ierr = MatPreallocated(mat);CHKERRQ(ierr); ierr = PetscLogEventBegin(MAT_CholeskyFactorSymbolic,mat,perm,0,0);CHKERRQ(ierr); ierr = (*mat->ops->choleskyfactorsymbolic)(mat,perm,info,fact);CHKERRQ(ierr); ierr = PetscLogEventEnd(MAT_CholeskyFactorSymbolic,mat,perm,0,0);CHKERRQ(ierr); ierr = PetscObjectStateIncrease((PetscObject)*fact);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatCholeskyFactorNumeric" /*@ MatCholeskyFactorNumeric - Performs numeric Cholesky factorization of a symmetric matrix. Call this routine after first calling MatCholeskyFactorSymbolic(). Collective on Mat Input Parameter: . mat - the initial matrix . info - options for factorization - fact - the symbolic factor of mat Output Parameter: . fact - the factored matrix Notes: Most users should employ the simplified KSP interface for linear solvers instead of working directly with matrix algebra routines such as this. See, e.g., KSPCreate(). Level: developer Concepts: matrices^Cholesky numeric factorization .seealso: MatCholeskyFactorSymbolic(), MatCholeskyFactor(), MatLUFactorNumeric() @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactorNumeric(Mat mat,MatFactorInfo *info,Mat *fact) { PetscErrorCode ierr; PetscFunctionBegin; PetscValidHeaderSpecific(mat,MAT_COOKIE,1); PetscValidType(mat,1); PetscValidPointer(fact,2); PetscValidHeaderSpecific(*fact,MAT_COOKIE,2); if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); if (!(*fact)->ops->choleskyfactornumeric) SETERRQ1(PETSC_ERR_SUP,"Mat type %s",mat->type_name); if (mat->rmap.N != (*fact)->rmap.N || mat->cmap.N != (*fact)->cmap.N) { SETERRQ4(PETSC_ERR_ARG_SIZ,"Mat mat,Mat *fact: global dim %D should = %D %D should = %D",mat->rmap.N,(*fact)->rmap.N,mat->cmap.N,(*fact)->cmap.N); } ierr = MatPreallocated(mat);CHKERRQ(ierr); ierr = PetscLogEventBegin(MAT_CholeskyFactorNumeric,mat,*fact,0,0);CHKERRQ(ierr); ierr = (*(*fact)->ops->choleskyfactornumeric)(mat,info,fact);CHKERRQ(ierr); ierr = PetscLogEventEnd(MAT_CholeskyFactorNumeric,mat,*fact,0,0);CHKERRQ(ierr); ierr = MatView_Private(*fact);CHKERRQ(ierr); ierr = PetscObjectStateIncrease((PetscObject)*fact);CHKERRQ(ierr); PetscFunctionReturn(0); } /* ----------------------------------------------------------------*/ #undef __FUNCT__ #define __FUNCT__ "MatSolve" /*@ MatSolve - Solves A x = b, given a factored matrix. Collective on Mat and Vec Input Parameters: + mat - the factored matrix - b - the right-hand-side vector Output Parameter: . x - the result vector Notes: The vectors b and x cannot be the same. I.e., one cannot call MatSolve(A,x,x). Notes: Most users should employ the simplified KSP interface for linear solvers instead of working directly with matrix algebra routines such as this. See, e.g., KSPCreate(). Level: developer Concepts: matrices^triangular solves .seealso: MatSolveAdd(), MatSolveTranspose(), MatSolveTransposeAdd() @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatSolve(Mat mat,Vec b,Vec x) { PetscErrorCode ierr; PetscFunctionBegin; PetscValidHeaderSpecific(mat,MAT_COOKIE,1); PetscValidType(mat,1); PetscValidHeaderSpecific(b,VEC_COOKIE,2); PetscValidHeaderSpecific(x,VEC_COOKIE,3); PetscCheckSameComm(mat,1,b,2); PetscCheckSameComm(mat,1,x,3); if (x == b) SETERRQ(PETSC_ERR_ARG_IDN,"x and b must be different vectors"); if (!mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Unfactored matrix"); if (mat->cmap.N != x->map.N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec x: global dim %D %D",mat->cmap.N,x->map.N); if (mat->rmap.N != b->map.N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec b: global dim %D %D",mat->rmap.N,b->map.N); if (mat->rmap.n != b->map.n) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec b: local dim %D %D",mat->rmap.n,b->map.n); if (!mat->rmap.N && !mat->cmap.N) PetscFunctionReturn(0); if (!mat->ops->solve) SETERRQ1(PETSC_ERR_SUP,"Mat type %s",mat->type_name); ierr = MatPreallocated(mat);CHKERRQ(ierr); ierr = PetscLogEventBegin(MAT_Solve,mat,b,x,0);CHKERRQ(ierr); ierr = (*mat->ops->solve)(mat,b,x);CHKERRQ(ierr); ierr = PetscLogEventEnd(MAT_Solve,mat,b,x,0);CHKERRQ(ierr); ierr = PetscObjectStateIncrease((PetscObject)x);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatMatSolve" /*@ MatMatSolve - Solves A X = B, given a factored matrix. Collective on Mat Input Parameters: + mat - the factored matrix - b - the right-hand-side vector Output Parameter: . x - the result vector Notes: The vectors b and x cannot be the same. I.e., one cannot call MatMatSolve(A,x,x). Notes: Most users should employ the simplified KSP interface for linear solvers instead of working directly with matrix algebra routines such as this. See, e.g., KSPCreate(). Level: developer Concepts: matrices^triangular solves .seealso: MatMatSolveAdd(), MatMatSolveTranspose(), MatMatSolveTransposeAdd() @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatMatSolve(Mat A,Mat B,Mat X) { PetscErrorCode ierr; PetscFunctionBegin; PetscValidHeaderSpecific(A,MAT_COOKIE,1); PetscValidType(A,1); PetscValidHeaderSpecific(B,MAT_COOKIE,2); PetscValidHeaderSpecific(X,MAT_COOKIE,3); PetscCheckSameComm(A,1,B,2); PetscCheckSameComm(A,1,X,3); if (X == B) SETERRQ(PETSC_ERR_ARG_IDN,"X and B must be different matrices"); if (!A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Unfactored matrix"); if (A->cmap.N != X->rmap.N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat A,Mat X: global dim %D %D",A->cmap.N,X->rmap.N); if (A->rmap.N != B->rmap.N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat A,Mat B: global dim %D %D",A->rmap.N,B->rmap.N); if (A->rmap.n != B->rmap.n) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat A,Mat B: local dim %D %D",A->rmap.n,B->rmap.n); if (!A->rmap.N && !A->cmap.N) PetscFunctionReturn(0); if (!A->ops->matsolve) SETERRQ1(PETSC_ERR_SUP,"Mat type %s",A->type_name); ierr = MatPreallocated(A);CHKERRQ(ierr); ierr = PetscLogEventBegin(MAT_MatSolve,A,B,X,0);CHKERRQ(ierr); ierr = (*A->ops->matsolve)(A,B,X);CHKERRQ(ierr); ierr = PetscLogEventEnd(MAT_MatSolve,A,B,X,0);CHKERRQ(ierr); ierr = PetscObjectStateIncrease((PetscObject)X);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatForwardSolve" /* @ MatForwardSolve - Solves L x = b, given a factored matrix, A = LU, or U^T*D^(1/2) x = b, given a factored symmetric matrix, A = U^T*D*U, Collective on Mat and Vec Input Parameters: + mat - the factored matrix - b - the right-hand-side vector Output Parameter: . x - the result vector Notes: MatSolve() should be used for most applications, as it performs a forward solve followed by a backward solve. The vectors b and x cannot be the same, i.e., one cannot call MatForwardSolve(A,x,x). For matrix in seqsbaij format with block size larger than 1, the diagonal blocks are not implemented as D = D^(1/2) * D^(1/2) yet. MatForwardSolve() solves U^T*D y = b, and MatBackwardSolve() solves U x = y. Thus they do not provide a symmetric preconditioner. Most users should employ the simplified KSP interface for linear solvers instead of working directly with matrix algebra routines such as this. See, e.g., KSPCreate(). Level: developer Concepts: matrices^forward solves .seealso: MatSolve(), MatBackwardSolve() @ */ PetscErrorCode PETSCMAT_DLLEXPORT MatForwardSolve(Mat mat,Vec b,Vec x) { PetscErrorCode ierr; PetscFunctionBegin; PetscValidHeaderSpecific(mat,MAT_COOKIE,1); PetscValidType(mat,1); PetscValidHeaderSpecific(b,VEC_COOKIE,2); PetscValidHeaderSpecific(x,VEC_COOKIE,3); PetscCheckSameComm(mat,1,b,2); PetscCheckSameComm(mat,1,x,3); if (x == b) SETERRQ(PETSC_ERR_ARG_IDN,"x and b must be different vectors"); if (!mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Unfactored matrix"); if (!mat->ops->forwardsolve) SETERRQ1(PETSC_ERR_SUP,"Mat type %s",mat->type_name); if (mat->cmap.N != x->map.N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec x: global dim %D %D",mat->cmap.N,x->map.N); if (mat->rmap.N != b->map.N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec b: global dim %D %D",mat->rmap.N,b->map.N); if (mat->rmap.n != b->map.n) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec b: local dim %D %D",mat->rmap.n,b->map.n); ierr = MatPreallocated(mat);CHKERRQ(ierr); ierr = PetscLogEventBegin(MAT_ForwardSolve,mat,b,x,0);CHKERRQ(ierr); ierr = (*mat->ops->forwardsolve)(mat,b,x);CHKERRQ(ierr); ierr = PetscLogEventEnd(MAT_ForwardSolve,mat,b,x,0);CHKERRQ(ierr); ierr = PetscObjectStateIncrease((PetscObject)x);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatBackwardSolve" /* @ MatBackwardSolve - Solves U x = b, given a factored matrix, A = LU. D^(1/2) U x = b, given a factored symmetric matrix, A = U^T*D*U, Collective on Mat and Vec Input Parameters: + mat - the factored matrix - b - the right-hand-side vector Output Parameter: . x - the result vector Notes: MatSolve() should be used for most applications, as it performs a forward solve followed by a backward solve. The vectors b and x cannot be the same. I.e., one cannot call MatBackwardSolve(A,x,x). For matrix in seqsbaij format with block size larger than 1, the diagonal blocks are not implemented as D = D^(1/2) * D^(1/2) yet. MatForwardSolve() solves U^T*D y = b, and MatBackwardSolve() solves U x = y. Thus they do not provide a symmetric preconditioner. Most users should employ the simplified KSP interface for linear solvers instead of working directly with matrix algebra routines such as this. See, e.g., KSPCreate(). Level: developer Concepts: matrices^backward solves .seealso: MatSolve(), MatForwardSolve() @ */ PetscErrorCode PETSCMAT_DLLEXPORT MatBackwardSolve(Mat mat,Vec b,Vec x) { PetscErrorCode ierr; PetscFunctionBegin; PetscValidHeaderSpecific(mat,MAT_COOKIE,1); PetscValidType(mat,1); PetscValidHeaderSpecific(b,VEC_COOKIE,2); PetscValidHeaderSpecific(x,VEC_COOKIE,3); PetscCheckSameComm(mat,1,b,2); PetscCheckSameComm(mat,1,x,3); if (x == b) SETERRQ(PETSC_ERR_ARG_IDN,"x and b must be different vectors"); if (!mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Unfactored matrix"); if (!mat->ops->backwardsolve) SETERRQ1(PETSC_ERR_SUP,"Mat type %s",mat->type_name); if (mat->cmap.N != x->map.N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec x: global dim %D %D",mat->cmap.N,x->map.N); if (mat->rmap.N != b->map.N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec b: global dim %D %D",mat->rmap.N,b->map.N); if (mat->rmap.n != b->map.n) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec b: local dim %D %D",mat->rmap.n,b->map.n); ierr = MatPreallocated(mat);CHKERRQ(ierr); ierr = PetscLogEventBegin(MAT_BackwardSolve,mat,b,x,0);CHKERRQ(ierr); ierr = (*mat->ops->backwardsolve)(mat,b,x);CHKERRQ(ierr); ierr = PetscLogEventEnd(MAT_BackwardSolve,mat,b,x,0);CHKERRQ(ierr); ierr = PetscObjectStateIncrease((PetscObject)x);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatSolveAdd" /*@ MatSolveAdd - Computes x = y + inv(A)*b, given a factored matrix. Collective on Mat and Vec Input Parameters: + mat - the factored matrix . b - the right-hand-side vector - y - the vector to be added to Output Parameter: . x - the result vector Notes: The vectors b and x cannot be the same. I.e., one cannot call MatSolveAdd(A,x,y,x). Most users should employ the simplified KSP interface for linear solvers instead of working directly with matrix algebra routines such as this. See, e.g., KSPCreate(). Level: developer Concepts: matrices^triangular solves .seealso: MatSolve(), MatSolveTranspose(), MatSolveTransposeAdd() @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatSolveAdd(Mat mat,Vec b,Vec y,Vec x) { PetscScalar one = 1.0; Vec tmp; PetscErrorCode ierr; PetscFunctionBegin; PetscValidHeaderSpecific(mat,MAT_COOKIE,1); PetscValidType(mat,1); PetscValidHeaderSpecific(y,VEC_COOKIE,2); PetscValidHeaderSpecific(b,VEC_COOKIE,3); PetscValidHeaderSpecific(x,VEC_COOKIE,4); PetscCheckSameComm(mat,1,b,2); PetscCheckSameComm(mat,1,y,2); PetscCheckSameComm(mat,1,x,3); if (x == b) SETERRQ(PETSC_ERR_ARG_IDN,"x and b must be different vectors"); if (!mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Unfactored matrix"); if (mat->cmap.N != x->map.N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec x: global dim %D %D",mat->cmap.N,x->map.N); if (mat->rmap.N != b->map.N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec b: global dim %D %D",mat->rmap.N,b->map.N); if (mat->rmap.N != y->map.N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec y: global dim %D %D",mat->rmap.N,y->map.N); if (mat->rmap.n != b->map.n) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec b: local dim %D %D",mat->rmap.n,b->map.n); if (x->map.n != y->map.n) SETERRQ2(PETSC_ERR_ARG_SIZ,"Vec x,Vec y: local dim %D %D",x->map.n,y->map.n); ierr = MatPreallocated(mat);CHKERRQ(ierr); ierr = PetscLogEventBegin(MAT_SolveAdd,mat,b,x,y);CHKERRQ(ierr); if (mat->ops->solveadd) { ierr = (*mat->ops->solveadd)(mat,b,y,x);CHKERRQ(ierr); } else { /* do the solve then the add manually */ if (x != y) { ierr = MatSolve(mat,b,x);CHKERRQ(ierr); ierr = VecAXPY(x,one,y);CHKERRQ(ierr); } else { ierr = VecDuplicate(x,&tmp);CHKERRQ(ierr); ierr = PetscLogObjectParent(mat,tmp);CHKERRQ(ierr); ierr = VecCopy(x,tmp);CHKERRQ(ierr); ierr = MatSolve(mat,b,x);CHKERRQ(ierr); ierr = VecAXPY(x,one,tmp);CHKERRQ(ierr); ierr = VecDestroy(tmp);CHKERRQ(ierr); } } ierr = PetscLogEventEnd(MAT_SolveAdd,mat,b,x,y);CHKERRQ(ierr); ierr = PetscObjectStateIncrease((PetscObject)x);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatSolveTranspose" /*@ MatSolveTranspose - Solves A' x = b, given a factored matrix. Collective on Mat and Vec Input Parameters: + mat - the factored matrix - b - the right-hand-side vector Output Parameter: . x - the result vector Notes: The vectors b and x cannot be the same. I.e., one cannot call MatSolveTranspose(A,x,x). Most users should employ the simplified KSP interface for linear solvers instead of working directly with matrix algebra routines such as this. See, e.g., KSPCreate(). Level: developer Concepts: matrices^triangular solves .seealso: MatSolve(), MatSolveAdd(), MatSolveTransposeAdd() @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatSolveTranspose(Mat mat,Vec b,Vec x) { PetscErrorCode ierr; PetscFunctionBegin; PetscValidHeaderSpecific(mat,MAT_COOKIE,1); PetscValidType(mat,1); PetscValidHeaderSpecific(b,VEC_COOKIE,2); PetscValidHeaderSpecific(x,VEC_COOKIE,3); PetscCheckSameComm(mat,1,b,2); PetscCheckSameComm(mat,1,x,3); if (!mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Unfactored matrix"); if (x == b) SETERRQ(PETSC_ERR_ARG_IDN,"x and b must be different vectors"); if (!mat->ops->solvetranspose) SETERRQ1(PETSC_ERR_SUP,"Matrix type %s",mat->type_name); if (mat->rmap.N != x->map.N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec x: global dim %D %D",mat->rmap.N,x->map.N); if (mat->cmap.N != b->map.N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec b: global dim %D %D",mat->cmap.N,b->map.N); ierr = MatPreallocated(mat);CHKERRQ(ierr); ierr = PetscLogEventBegin(MAT_SolveTranspose,mat,b,x,0);CHKERRQ(ierr); ierr = (*mat->ops->solvetranspose)(mat,b,x);CHKERRQ(ierr); ierr = PetscLogEventEnd(MAT_SolveTranspose,mat,b,x,0);CHKERRQ(ierr); ierr = PetscObjectStateIncrease((PetscObject)x);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatSolveTransposeAdd" /*@ MatSolveTransposeAdd - Computes x = y + inv(Transpose(A)) b, given a factored matrix. Collective on Mat and Vec Input Parameters: + mat - the factored matrix . b - the right-hand-side vector - y - the vector to be added to Output Parameter: . x - the result vector Notes: The vectors b and x cannot be the same. I.e., one cannot call MatSolveTransposeAdd(A,x,y,x). Most users should employ the simplified KSP interface for linear solvers instead of working directly with matrix algebra routines such as this. See, e.g., KSPCreate(). Level: developer Concepts: matrices^triangular solves .seealso: MatSolve(), MatSolveAdd(), MatSolveTranspose() @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatSolveTransposeAdd(Mat mat,Vec b,Vec y,Vec x) { PetscScalar one = 1.0; PetscErrorCode ierr; Vec tmp; PetscFunctionBegin; PetscValidHeaderSpecific(mat,MAT_COOKIE,1); PetscValidType(mat,1); PetscValidHeaderSpecific(y,VEC_COOKIE,2); PetscValidHeaderSpecific(b,VEC_COOKIE,3); PetscValidHeaderSpecific(x,VEC_COOKIE,4); PetscCheckSameComm(mat,1,b,2); PetscCheckSameComm(mat,1,y,3); PetscCheckSameComm(mat,1,x,4); if (x == b) SETERRQ(PETSC_ERR_ARG_IDN,"x and b must be different vectors"); if (!mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Unfactored matrix"); if (mat->rmap.N != x->map.N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec x: global dim %D %D",mat->rmap.N,x->map.N); if (mat->cmap.N != b->map.N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec b: global dim %D %D",mat->cmap.N,b->map.N); if (mat->cmap.N != y->map.N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec y: global dim %D %D",mat->cmap.N,y->map.N); if (x->map.n != y->map.n) SETERRQ2(PETSC_ERR_ARG_SIZ,"Vec x,Vec y: local dim %D %D",x->map.n,y->map.n); ierr = MatPreallocated(mat);CHKERRQ(ierr); ierr = PetscLogEventBegin(MAT_SolveTransposeAdd,mat,b,x,y);CHKERRQ(ierr); if (mat->ops->solvetransposeadd) { ierr = (*mat->ops->solvetransposeadd)(mat,b,y,x);CHKERRQ(ierr); } else { /* do the solve then the add manually */ if (x != y) { ierr = MatSolveTranspose(mat,b,x);CHKERRQ(ierr); ierr = VecAXPY(x,one,y);CHKERRQ(ierr); } else { ierr = VecDuplicate(x,&tmp);CHKERRQ(ierr); ierr = PetscLogObjectParent(mat,tmp);CHKERRQ(ierr); ierr = VecCopy(x,tmp);CHKERRQ(ierr); ierr = MatSolveTranspose(mat,b,x);CHKERRQ(ierr); ierr = VecAXPY(x,one,tmp);CHKERRQ(ierr); ierr = VecDestroy(tmp);CHKERRQ(ierr); } } ierr = PetscLogEventEnd(MAT_SolveTransposeAdd,mat,b,x,y);CHKERRQ(ierr); ierr = PetscObjectStateIncrease((PetscObject)x);CHKERRQ(ierr); PetscFunctionReturn(0); } /* ----------------------------------------------------------------*/ #undef __FUNCT__ #define __FUNCT__ "MatRelax" /*@ MatRelax - Computes relaxation (SOR, Gauss-Seidel) sweeps. Collective on Mat and Vec Input Parameters: + mat - the matrix . b - the right hand side . omega - the relaxation factor . flag - flag indicating the type of SOR (see below) . shift - diagonal shift - its - the number of iterations - lits - the number of local iterations Output Parameters: . x - the solution (can contain an initial guess) SOR Flags: . SOR_FORWARD_SWEEP - forward SOR . SOR_BACKWARD_SWEEP - backward SOR . SOR_SYMMETRIC_SWEEP - SSOR (symmetric SOR) . SOR_LOCAL_FORWARD_SWEEP - local forward SOR . SOR_LOCAL_BACKWARD_SWEEP - local forward SOR . SOR_LOCAL_SYMMETRIC_SWEEP - local SSOR . SOR_APPLY_UPPER, SOR_APPLY_LOWER - applies upper/lower triangular part of matrix to vector (with omega) . SOR_ZERO_INITIAL_GUESS - zero initial guess Notes: SOR_LOCAL_FORWARD_SWEEP, SOR_LOCAL_BACKWARD_SWEEP, and SOR_LOCAL_SYMMETRIC_SWEEP perform separate independent smoothings on each processor. Application programmers will not generally use MatRelax() directly, but instead will employ the KSP/PC interface. Notes for Advanced Users: The flags are implemented as bitwise inclusive or operations. For example, use (SOR_ZERO_INITIAL_GUESS | SOR_SYMMETRIC_SWEEP) to specify a zero initial guess for SSOR. Most users should employ the simplified KSP interface for linear solvers instead of working directly with matrix algebra routines such as this. See, e.g., KSPCreate(). See also, MatPBRelax(). This routine will automatically call the point block version if the point version is not available. Level: developer Concepts: matrices^relaxation Concepts: matrices^SOR Concepts: matrices^Gauss-Seidel @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatRelax(Mat mat,Vec b,PetscReal omega,MatSORType flag,PetscReal shift,PetscInt its,PetscInt lits,Vec x) { PetscErrorCode ierr; PetscFunctionBegin; PetscValidHeaderSpecific(mat,MAT_COOKIE,1); PetscValidType(mat,1); PetscValidHeaderSpecific(b,VEC_COOKIE,2); PetscValidHeaderSpecific(x,VEC_COOKIE,8); PetscCheckSameComm(mat,1,b,2); PetscCheckSameComm(mat,1,x,8); if (!mat->ops->relax && !mat->ops->pbrelax) SETERRQ1(PETSC_ERR_SUP,"Mat type %s",mat->type_name); if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); if (mat->cmap.N != x->map.N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec x: global dim %D %D",mat->cmap.N,x->map.N); if (mat->rmap.N != b->map.N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec b: global dim %D %D",mat->rmap.N,b->map.N); if (mat->rmap.n != b->map.n) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec b: local dim %D %D",mat->rmap.n,b->map.n); ierr = MatPreallocated(mat);CHKERRQ(ierr); ierr = PetscLogEventBegin(MAT_Relax,mat,b,x,0);CHKERRQ(ierr); if (mat->ops->relax) { ierr =(*mat->ops->relax)(mat,b,omega,flag,shift,its,lits,x);CHKERRQ(ierr); } else { ierr =(*mat->ops->pbrelax)(mat,b,omega,flag,shift,its,lits,x);CHKERRQ(ierr); } ierr = PetscLogEventEnd(MAT_Relax,mat,b,x,0);CHKERRQ(ierr); ierr = PetscObjectStateIncrease((PetscObject)x);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatPBRelax" /*@ MatPBRelax - Computes relaxation (SOR, Gauss-Seidel) sweeps. Collective on Mat and Vec See MatRelax() for usage For multi-component PDEs where the Jacobian is stored in a point block format (with the PETSc BAIJ matrix formats) the relaxation is done one point block at a time. That is, the small (for example, 4 by 4) blocks along the diagonal are solved simultaneously (that is a 4 by 4 linear solve is done) to update all the values at a point. Level: developer @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatPBRelax(Mat mat,Vec b,PetscReal omega,MatSORType flag,PetscReal shift,PetscInt its,PetscInt lits,Vec x) { PetscErrorCode ierr; PetscFunctionBegin; PetscValidHeaderSpecific(mat,MAT_COOKIE,1); PetscValidType(mat,1); PetscValidHeaderSpecific(b,VEC_COOKIE,2); PetscValidHeaderSpecific(x,VEC_COOKIE,8); PetscCheckSameComm(mat,1,b,2); PetscCheckSameComm(mat,1,x,8); if (!mat->ops->pbrelax) SETERRQ1(PETSC_ERR_SUP,"Mat type %s",mat->type_name); if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); if (mat->cmap.N != x->map.N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec x: global dim %D %D",mat->cmap.N,x->map.N); if (mat->rmap.N != b->map.N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec b: global dim %D %D",mat->rmap.N,b->map.N); if (mat->rmap.n != b->map.n) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec b: local dim %D %D",mat->rmap.n,b->map.n); ierr = MatPreallocated(mat);CHKERRQ(ierr); ierr = PetscLogEventBegin(MAT_Relax,mat,b,x,0);CHKERRQ(ierr); ierr =(*mat->ops->pbrelax)(mat,b,omega,flag,shift,its,lits,x);CHKERRQ(ierr); ierr = PetscLogEventEnd(MAT_Relax,mat,b,x,0);CHKERRQ(ierr); ierr = PetscObjectStateIncrease((PetscObject)x);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatCopy_Basic" /* Default matrix copy routine. */ PetscErrorCode MatCopy_Basic(Mat A,Mat B,MatStructure str) { PetscErrorCode ierr; PetscInt i,rstart,rend,nz; const PetscInt *cwork; const PetscScalar *vwork; PetscFunctionBegin; if (B->assembled) { ierr = MatZeroEntries(B);CHKERRQ(ierr); } ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); for (i=rstart; iassembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); if (A->rmap.N != B->rmap.N || A->cmap.N != B->cmap.N) SETERRQ4(PETSC_ERR_ARG_SIZ,"Mat A,Mat B: global dim (%D,%D) (%D,%D)",A->rmap.N,B->rmap.N,A->cmap.N,B->cmap.N); ierr = MatPreallocated(A);CHKERRQ(ierr); ierr = PetscLogEventBegin(MAT_Copy,A,B,0,0);CHKERRQ(ierr); if (A->ops->copy) { ierr = (*A->ops->copy)(A,B,str);CHKERRQ(ierr); } else { /* generic conversion */ ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr); } if (A->mapping) { if (B->mapping) {ierr = ISLocalToGlobalMappingDestroy(B->mapping);CHKERRQ(ierr);B->mapping = 0;} ierr = MatSetLocalToGlobalMapping(B,A->mapping);CHKERRQ(ierr); } if (A->bmapping) { if (B->bmapping) {ierr = ISLocalToGlobalMappingDestroy(B->bmapping);CHKERRQ(ierr);B->bmapping = 0;} ierr = MatSetLocalToGlobalMappingBlock(B,A->mapping);CHKERRQ(ierr); } ierr = PetscLogEventEnd(MAT_Copy,A,B,0,0);CHKERRQ(ierr); ierr = PetscObjectStateIncrease((PetscObject)B);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatConvert" /*@C MatConvert - Converts a matrix to another matrix, either of the same or different type. Collective on Mat Input Parameters: + mat - the matrix . newtype - new matrix type. Use MATSAME to create a new matrix of the same type as the original matrix. - reuse - denotes if the destination matrix is to be created or reused. Currently MAT_REUSE_MATRIX is only supported for inplace conversion, otherwise use MAT_INITIAL_MATRIX. Output Parameter: . M - pointer to place new matrix Notes: MatConvert() first creates a new matrix and then copies the data from the first matrix. A related routine is MatCopy(), which copies the matrix entries of one matrix to another already existing matrix context. Level: intermediate Concepts: matrices^converting between storage formats .seealso: MatCopy(), MatDuplicate() @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatConvert(Mat mat, MatType newtype,MatReuse reuse,Mat *M) { PetscErrorCode ierr; PetscTruth sametype,issame,flg; char convname[256],mtype[256]; Mat B; PetscFunctionBegin; PetscValidHeaderSpecific(mat,MAT_COOKIE,1); PetscValidType(mat,1); PetscValidPointer(M,3); if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); ierr = MatPreallocated(mat);CHKERRQ(ierr); ierr = PetscOptionsGetString(PETSC_NULL,"-matconvert_type",mtype,256,&flg);CHKERRQ(ierr); if (flg) { newtype = mtype; } ierr = PetscTypeCompare((PetscObject)mat,newtype,&sametype);CHKERRQ(ierr); ierr = PetscStrcmp(newtype,"same",&issame);CHKERRQ(ierr); if ((reuse==MAT_REUSE_MATRIX) && (mat != *M)) { SETERRQ(PETSC_ERR_SUP,"MAT_REUSE_MATRIX only supported for in-place conversion currently"); } if ((sametype || issame) && (reuse==MAT_INITIAL_MATRIX) && mat->ops->duplicate) { ierr = (*mat->ops->duplicate)(mat,MAT_COPY_VALUES,M);CHKERRQ(ierr); } else { PetscErrorCode (*conv)(Mat, MatType,MatReuse,Mat*)=PETSC_NULL; const char *prefix[3] = {"seq","mpi",""}; PetscInt i; /* Order of precedence: 1) See if a specialized converter is known to the current matrix. 2) See if a specialized converter is known to the desired matrix class. 3) See if a good general converter is registered for the desired class (as of 6/27/03 only MATMPIADJ falls into this category). 4) See if a good general converter is known for the current matrix. 5) Use a really basic converter. */ /* 1) See if a specialized converter is known to the current matrix and the desired class */ for (i=0; i<3; i++) { ierr = PetscStrcpy(convname,"MatConvert_");CHKERRQ(ierr); ierr = PetscStrcat(convname,mat->type_name);CHKERRQ(ierr); ierr = PetscStrcat(convname,"_");CHKERRQ(ierr); ierr = PetscStrcat(convname,prefix[i]);CHKERRQ(ierr); ierr = PetscStrcat(convname,newtype);CHKERRQ(ierr); ierr = PetscStrcat(convname,"_C");CHKERRQ(ierr); ierr = PetscObjectQueryFunction((PetscObject)mat,convname,(void (**)(void))&conv);CHKERRQ(ierr); if (conv) goto foundconv; } /* 2) See if a specialized converter is known to the desired matrix class. */ ierr = MatCreate(mat->comm,&B);CHKERRQ(ierr); ierr = MatSetSizes(B,mat->rmap.n,mat->cmap.n,mat->rmap.N,mat->cmap.N);CHKERRQ(ierr); ierr = MatSetType(B,newtype);CHKERRQ(ierr); for (i=0; i<3; i++) { ierr = PetscStrcpy(convname,"MatConvert_");CHKERRQ(ierr); ierr = PetscStrcat(convname,mat->type_name);CHKERRQ(ierr); ierr = PetscStrcat(convname,"_");CHKERRQ(ierr); ierr = PetscStrcat(convname,prefix[i]);CHKERRQ(ierr); ierr = PetscStrcat(convname,newtype);CHKERRQ(ierr); ierr = PetscStrcat(convname,"_C");CHKERRQ(ierr); ierr = PetscObjectQueryFunction((PetscObject)B,convname,(void (**)(void))&conv);CHKERRQ(ierr); if (conv) { ierr = MatDestroy(B);CHKERRQ(ierr); goto foundconv; } } /* 3) See if a good general converter is registered for the desired class */ conv = B->ops->convertfrom; ierr = MatDestroy(B);CHKERRQ(ierr); if (conv) goto foundconv; /* 4) See if a good general converter is known for the current matrix */ if (mat->ops->convert) { conv = mat->ops->convert; } if (conv) goto foundconv; /* 5) Use a really basic converter. */ conv = MatConvert_Basic; foundconv: ierr = PetscLogEventBegin(MAT_Convert,mat,0,0,0);CHKERRQ(ierr); ierr = (*conv)(mat,newtype,reuse,M);CHKERRQ(ierr); ierr = PetscLogEventEnd(MAT_Convert,mat,0,0,0);CHKERRQ(ierr); } B = *M; ierr = PetscObjectStateIncrease((PetscObject)B);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatDuplicate" /*@ MatDuplicate - Duplicates a matrix including the non-zero structure. Collective on Mat Input Parameters: + mat - the matrix - op - either MAT_DO_NOT_COPY_VALUES or MAT_COPY_VALUES, cause it to copy nonzero values as well or not Output Parameter: . M - pointer to place new matrix Level: intermediate Concepts: matrices^duplicating .seealso: MatCopy(), MatConvert() @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatDuplicate(Mat mat,MatDuplicateOption op,Mat *M) { PetscErrorCode ierr; Mat B; PetscFunctionBegin; PetscValidHeaderSpecific(mat,MAT_COOKIE,1); PetscValidType(mat,1); PetscValidPointer(M,3); if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); ierr = MatPreallocated(mat);CHKERRQ(ierr); *M = 0; if (!mat->ops->duplicate) { SETERRQ(PETSC_ERR_SUP,"Not written for this matrix type"); } ierr = PetscLogEventBegin(MAT_Convert,mat,0,0,0);CHKERRQ(ierr); ierr = (*mat->ops->duplicate)(mat,op,M);CHKERRQ(ierr); B = *M; if (mat->mapping) { ierr = MatSetLocalToGlobalMapping(B,mat->mapping);CHKERRQ(ierr); } if (mat->bmapping) { ierr = MatSetLocalToGlobalMappingBlock(B,mat->bmapping);CHKERRQ(ierr); } ierr = PetscMapCopy(mat->comm,&mat->rmap,&B->rmap);CHKERRQ(ierr); ierr = PetscMapCopy(mat->comm,&mat->cmap,&B->cmap);CHKERRQ(ierr); ierr = PetscLogEventEnd(MAT_Convert,mat,0,0,0);CHKERRQ(ierr); ierr = PetscObjectStateIncrease((PetscObject)B);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatGetDiagonal" /*@ MatGetDiagonal - Gets the diagonal of a matrix. Collective on Mat and Vec Input Parameters: + mat - the matrix - v - the vector for storing the diagonal Output Parameter: . v - the diagonal of the matrix Notes: The result of this call are the same as if one converted the matrix to dense format and found the minimum value in each row (i.e. the implicit zeros are counted as zeros). Level: intermediate Concepts: matrices^accessing diagonals .seealso: MatGetRow(), MatGetSubMatrices(), MatGetSubmatrix(), MatGetRowMaxAbs() @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatGetDiagonal(Mat mat,Vec v) { PetscErrorCode ierr; PetscFunctionBegin; PetscValidHeaderSpecific(mat,MAT_COOKIE,1); PetscValidType(mat,1); PetscValidHeaderSpecific(v,VEC_COOKIE,2); if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); if (!mat->ops->getdiagonal) SETERRQ1(PETSC_ERR_SUP,"Mat type %s",mat->type_name); ierr = MatPreallocated(mat);CHKERRQ(ierr); ierr = (*mat->ops->getdiagonal)(mat,v);CHKERRQ(ierr); ierr = PetscObjectStateIncrease((PetscObject)v);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatGetRowMin" /*@ MatGetRowMin - Gets the minimum value (of the real part) of each row of the matrix Collective on Mat and Vec Input Parameters: . mat - the matrix Output Parameter: + v - the vector for storing the maximums - idx - the indices of the column found for each row (optional) Level: intermediate Notes: The result of this call are the same as if one converted the matrix to dense format and found the minimum value in each row (i.e. the implicit zeros are counted as zeros). This code is only implemented for a couple of matrix formats. Concepts: matrices^getting row maximums .seealso: MatGetDiagonal(), MatGetSubMatrices(), MatGetSubmatrix(), MatGetRowMaxAbs(), MatGetRowMax() @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMin(Mat mat,Vec v,PetscInt idx[]) { PetscErrorCode ierr; PetscFunctionBegin; PetscValidHeaderSpecific(mat,MAT_COOKIE,1); PetscValidType(mat,1); PetscValidHeaderSpecific(v,VEC_COOKIE,2); if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); if (!mat->ops->getrowmax) SETERRQ1(PETSC_ERR_SUP,"Mat type %s",mat->type_name); ierr = MatPreallocated(mat);CHKERRQ(ierr); ierr = (*mat->ops->getrowmin)(mat,v,idx);CHKERRQ(ierr); ierr = PetscObjectStateIncrease((PetscObject)v);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatGetRowMax" /*@ MatGetRowMax - Gets the maximum value (of the real part) of each row of the matrix Collective on Mat and Vec Input Parameters: . mat - the matrix Output Parameter: + v - the vector for storing the maximums - idx - the indices of the column found for each row (optional) Level: intermediate Notes: The result of this call are the same as if one converted the matrix to dense format and found the minimum value in each row (i.e. the implicit zeros are counted as zeros). This code is only implemented for a couple of matrix formats. Concepts: matrices^getting row maximums .seealso: MatGetDiagonal(), MatGetSubMatrices(), MatGetSubmatrix(), MatGetRowMaxAbs(), MatGetRowMin() @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMax(Mat mat,Vec v,PetscInt idx[]) { PetscErrorCode ierr; PetscFunctionBegin; PetscValidHeaderSpecific(mat,MAT_COOKIE,1); PetscValidType(mat,1); PetscValidHeaderSpecific(v,VEC_COOKIE,2); if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); if (!mat->ops->getrowmax) SETERRQ1(PETSC_ERR_SUP,"Mat type %s",mat->type_name); ierr = MatPreallocated(mat);CHKERRQ(ierr); ierr = (*mat->ops->getrowmax)(mat,v,idx);CHKERRQ(ierr); ierr = PetscObjectStateIncrease((PetscObject)v);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatGetRowMaxAbs" /*@ MatGetRowMaxAbs - Gets the maximum value (in absolute value) of each row of the matrix Collective on Mat and Vec Input Parameters: . mat - the matrix Output Parameter: + v - the vector for storing the maximums - idx - the indices of the column found for each row (optional) Level: intermediate Notes: if a row is completely empty or has only 0.0 values then the idx[] value for that row is 0 (the first column). This code is only implemented for a couple of matrix formats. Concepts: matrices^getting row maximums .seealso: MatGetDiagonal(), MatGetSubMatrices(), MatGetSubmatrix(), MatGetRowMax(), MatGetRowMin() @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMaxAbs(Mat mat,Vec v,PetscInt idx[]) { PetscErrorCode ierr; PetscFunctionBegin; PetscValidHeaderSpecific(mat,MAT_COOKIE,1); PetscValidType(mat,1); PetscValidHeaderSpecific(v,VEC_COOKIE,2); if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); if (!mat->ops->getrowmaxabs) SETERRQ1(PETSC_ERR_SUP,"Mat type %s",mat->type_name); ierr = MatPreallocated(mat);CHKERRQ(ierr); ierr = (*mat->ops->getrowmaxabs)(mat,v,idx);CHKERRQ(ierr); ierr = PetscObjectStateIncrease((PetscObject)v);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatGetRowSum" /*@ MatGetRowSum - Gets the sum of each row of the matrix Collective on Mat and Vec Input Parameters: . mat - the matrix Output Parameter: . v - the vector for storing the maximums Level: intermediate Notes: This code is slow since it is not currently specialized for different formats Concepts: matrices^getting row sums .seealso: MatGetDiagonal(), MatGetSubMatrices(), MatGetSubmatrix(), MatGetRowMax(), MatGetRowMin() @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowSum(Mat mat, Vec v) { PetscInt start, end, row; PetscScalar *array; PetscErrorCode ierr; PetscFunctionBegin; PetscValidHeaderSpecific(mat,MAT_COOKIE,1); PetscValidType(mat,1); PetscValidHeaderSpecific(v,VEC_COOKIE,2); if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); ierr = MatPreallocated(mat);CHKERRQ(ierr); ierr = MatGetOwnershipRange(mat, &start, &end);CHKERRQ(ierr); ierr = VecGetArray(v, &array);CHKERRQ(ierr); for(row = start; row < end; ++row) { PetscInt ncols, col; const PetscInt *cols; const PetscScalar *vals; array[row - start] = 0.0; ierr = MatGetRow(mat, row, &ncols, &cols, &vals);CHKERRQ(ierr); for(col = 0; col < ncols; col++) { array[row - start] += vals[col]; } } ierr = VecRestoreArray(v, &array);CHKERRQ(ierr); ierr = PetscObjectStateIncrease((PetscObject) v);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatTranspose" /*@C MatTranspose - Computes an in-place or out-of-place transpose of a matrix. Collective on Mat Input Parameter: . mat - the matrix to transpose Output Parameters: . B - the transpose Notes: If you pass in PETSC_NULL for B an in-place transpose in mat will be done Level: intermediate Concepts: matrices^transposing .seealso: MatMultTranspose(), MatMultTransposeAdd(), MatIsTranspose() @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatTranspose(Mat mat,Mat *B) { PetscErrorCode ierr; PetscFunctionBegin; PetscValidHeaderSpecific(mat,MAT_COOKIE,1); PetscValidType(mat,1); if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); if (!mat->ops->transpose) SETERRQ1(PETSC_ERR_SUP,"Mat type %s",mat->type_name); ierr = MatPreallocated(mat);CHKERRQ(ierr); ierr = PetscLogEventBegin(MAT_Transpose,mat,0,0,0);CHKERRQ(ierr); ierr = (*mat->ops->transpose)(mat,B);CHKERRQ(ierr); ierr = PetscLogEventEnd(MAT_Transpose,mat,0,0,0);CHKERRQ(ierr); if (B) {ierr = PetscObjectStateIncrease((PetscObject)*B);CHKERRQ(ierr);} PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatIsTranspose" /*@C MatIsTranspose - Test whether a matrix is another one's transpose, or its own, in which case it tests symmetry. Collective on Mat Input Parameter: + A - the matrix to test - B - the matrix to test against, this can equal the first parameter Output Parameters: . flg - the result Notes: Only available for SeqAIJ/MPIAIJ matrices. The sequential algorithm has a running time of the order of the number of nonzeros; the parallel test involves parallel copies of the block-offdiagonal parts of the matrix. Level: intermediate Concepts: matrices^transposing, matrix^symmetry .seealso: MatTranspose(), MatIsSymmetric(), MatIsHermitian() @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatIsTranspose(Mat A,Mat B,PetscReal tol,PetscTruth *flg) { PetscErrorCode ierr,(*f)(Mat,Mat,PetscReal,PetscTruth*),(*g)(Mat,Mat,PetscReal,PetscTruth*); PetscFunctionBegin; PetscValidHeaderSpecific(A,MAT_COOKIE,1); PetscValidHeaderSpecific(B,MAT_COOKIE,2); PetscValidPointer(flg,3); ierr = PetscObjectQueryFunction((PetscObject)A,"MatIsTranspose_C",(void (**)(void))&f);CHKERRQ(ierr); ierr = PetscObjectQueryFunction((PetscObject)B,"MatIsTranspose_C",(void (**)(void))&g);CHKERRQ(ierr); if (f && g) { if (f==g) { ierr = (*f)(A,B,tol,flg);CHKERRQ(ierr); } else { SETERRQ(PETSC_ERR_ARG_NOTSAMETYPE,"Matrices do not have the same comparator for symmetry test"); } } PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatPermute" /*@C MatPermute - Creates a new matrix with rows and columns permuted from the original. Collective on Mat Input Parameters: + mat - the matrix to permute . row - row permutation, each processor supplies only the permutation for its rows - col - column permutation, each processor needs the entire column permutation, that is this is the same size as the total number of columns in the matrix Output Parameters: . B - the permuted matrix Level: advanced Concepts: matrices^permuting .seealso: MatGetOrdering() @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatPermute(Mat mat,IS row,IS col,Mat *B) { PetscErrorCode ierr; PetscFunctionBegin; PetscValidHeaderSpecific(mat,MAT_COOKIE,1); PetscValidType(mat,1); PetscValidHeaderSpecific(row,IS_COOKIE,2); PetscValidHeaderSpecific(col,IS_COOKIE,3); PetscValidPointer(B,4); if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); if (!mat->ops->permute) SETERRQ1(PETSC_ERR_SUP,"MatPermute not available for Mat type %s",mat->type_name); ierr = MatPreallocated(mat);CHKERRQ(ierr); ierr = (*mat->ops->permute)(mat,row,col,B);CHKERRQ(ierr); ierr = PetscObjectStateIncrease((PetscObject)*B);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatPermuteSparsify" /*@C MatPermuteSparsify - Creates a new matrix with rows and columns permuted from the original and sparsified to the prescribed tolerance. Collective on Mat Input Parameters: + A - The matrix to permute . band - The half-bandwidth of the sparsified matrix, or PETSC_DECIDE . frac - The half-bandwidth as a fraction of the total size, or 0.0 . tol - The drop tolerance . rowp - The row permutation - colp - The column permutation Output Parameter: . B - The permuted, sparsified matrix Level: advanced Note: The default behavior (band = PETSC_DECIDE and frac = 0.0) is to restrict the half-bandwidth of the resulting matrix to 5% of the total matrix size. .keywords: matrix, permute, sparsify .seealso: MatGetOrdering(), MatPermute() @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatPermuteSparsify(Mat A, PetscInt band, PetscReal frac, PetscReal tol, IS rowp, IS colp, Mat *B) { IS irowp, icolp; PetscInt *rows, *cols; PetscInt M, N, locRowStart, locRowEnd; PetscInt nz, newNz; const PetscInt *cwork; PetscInt *cnew; const PetscScalar *vwork; PetscScalar *vnew; PetscInt bw, issize; PetscInt row, locRow, newRow, col, newCol; PetscErrorCode ierr; PetscFunctionBegin; PetscValidHeaderSpecific(A, MAT_COOKIE,1); PetscValidHeaderSpecific(rowp, IS_COOKIE,5); PetscValidHeaderSpecific(colp, IS_COOKIE,6); PetscValidPointer(B,7); if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); if (!A->ops->permutesparsify) { ierr = MatGetSize(A, &M, &N);CHKERRQ(ierr); ierr = MatGetOwnershipRange(A, &locRowStart, &locRowEnd);CHKERRQ(ierr); ierr = ISGetSize(rowp, &issize);CHKERRQ(ierr); if (issize != M) SETERRQ2(PETSC_ERR_ARG_WRONG, "Wrong size %D for row permutation, should be %D", issize, M); ierr = ISGetSize(colp, &issize);CHKERRQ(ierr); if (issize != N) SETERRQ2(PETSC_ERR_ARG_WRONG, "Wrong size %D for column permutation, should be %D", issize, N); ierr = ISInvertPermutation(rowp, 0, &irowp);CHKERRQ(ierr); ierr = ISGetIndices(irowp, &rows);CHKERRQ(ierr); ierr = ISInvertPermutation(colp, 0, &icolp);CHKERRQ(ierr); ierr = ISGetIndices(icolp, &cols);CHKERRQ(ierr); ierr = PetscMalloc(N * sizeof(PetscInt), &cnew);CHKERRQ(ierr); ierr = PetscMalloc(N * sizeof(PetscScalar), &vnew);CHKERRQ(ierr); /* Setup bandwidth to include */ if (band == PETSC_DECIDE) { if (frac <= 0.0) bw = (PetscInt) (M * 0.05); else bw = (PetscInt) (M * frac); } else { if (band <= 0) SETERRQ(PETSC_ERR_ARG_WRONG, "Bandwidth must be a positive integer"); bw = band; } /* Put values into new matrix */ ierr = MatDuplicate(A, MAT_DO_NOT_COPY_VALUES, B);CHKERRQ(ierr); for(row = locRowStart, locRow = 0; row < locRowEnd; row++, locRow++) { ierr = MatGetRow(A, row, &nz, &cwork, &vwork);CHKERRQ(ierr); newRow = rows[locRow]+locRowStart; for(col = 0, newNz = 0; col < nz; col++) { newCol = cols[cwork[col]]; if ((newCol >= newRow - bw) && (newCol < newRow + bw) && (PetscAbsScalar(vwork[col]) >= tol)) { cnew[newNz] = newCol; vnew[newNz] = vwork[col]; newNz++; } } ierr = MatSetValues(*B, 1, &newRow, newNz, cnew, vnew, INSERT_VALUES);CHKERRQ(ierr); ierr = MatRestoreRow(A, row, &nz, &cwork, &vwork);CHKERRQ(ierr); } ierr = PetscFree(cnew);CHKERRQ(ierr); ierr = PetscFree(vnew);CHKERRQ(ierr); ierr = MatAssemblyBegin(*B, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); ierr = MatAssemblyEnd(*B, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); ierr = ISRestoreIndices(irowp, &rows);CHKERRQ(ierr); ierr = ISRestoreIndices(icolp, &cols);CHKERRQ(ierr); ierr = ISDestroy(irowp);CHKERRQ(ierr); ierr = ISDestroy(icolp);CHKERRQ(ierr); } else { ierr = (*A->ops->permutesparsify)(A, band, frac, tol, rowp, co