[petsc-users] Some issues about DMDA and MatSetValuesStencil
Matthew Knepley
knepley at gmail.com
Wed Aug 3 09:39:45 CDT 2022
On Wed, Aug 3, 2022 at 9:17 AM wangxq2020--- via petsc-users <
petsc-users at mcs.anl.gov> wrote:
> Hello!
> I am a beginner about petsc,and I'm studing the DMDA recently. I hava a
> exercise code below, i run it with 4 processes, why is the size of A matrix
> is 15×15?
>
Are you sure it is not 16x16? The grid has 4 vertices in both the x- and
y-directions, so 16 vertices. You have 1 dof per vertex. This means that a
function
over this grid has 16 variables, so vectors have length 16. The Jacobian
relates the (linearized) change in the residual for a change in each
variables,
so it is a matrix of size 16 x 16.
Thanks,
Matt
> The user manual said A is the Jobian matix, is there any detailed material
> that explain the process? Another question is the setting values for A. I
> don't understand why a loop using the grid distribution information of each
> process can assign a value to A. The MatStencil mentioned in the user
> manual is a Data structure (C struct) for storing information about a
> single row or column of a matrix, why MatStencil contains i and j at the
> same time, and the row only needs i. Very much looking forward to
> reply,thanks.
>
> DM grid;
> ierr = DMDACreate2d // IN:
> ( comm, // collective on this communicator
> DM_BOUNDARY_NONE,DM_BOUNDARY_NONE, // no periodicity and such
> DMDA_STENCIL_STAR, // no cross connections
> 4,4, // global size 100x100; can be changed with
> options
> PETSC_DECIDE,PETSC_DECIDE, // processors in each direction
> 1, // degree of freedom per node
> 1, // stencil width
> NULL,NULL, // arrays of local sizes in each direction
> &grid // OUT: resulting object
> ); CHKERRQ(ierr);
> ierr = DMSetUp(grid); CHKERRQ(ierr);
>
> Mat A;
> ierr = DMCreateMatrix(grid,&A); CHKERRQ(ierr);
> DMDALocalInfo info;
> ierr = DMDAGetLocalInfo(grid,&info);CHKERRQ(ierr);
> for (int j=info.ys; j<info.ys+info.ym; j++) {
> for (int i=info.xs; i<info.xs+info.xm; i++) {
> MatStencil row = {0},col[5] = {{0}};
> PetscScalar v[5];
> PetscInt ncols = 0;
> row.j = j; row.i = i;
> printf("procno: %d;j:%d; i:%d\n",procno,j,i);
> /**** local connection: diagonal element ****/
> col[ncols].j = j; col[ncols].i = i; v[ncols++] = 4.;
> printf("procno: %d, j:%d. col[%d].j = %d, col[%d].i = %d,
> v[%d]=4\n",procno,j,ncols-1,j,ncols-1,i,ncols);
> /* boundaries: top row */
> if (i>0) {
> col[ncols].j = j; col[ncols].i = i-1; v[ncols++] = -1.;
> printf("procno: %d, j:%d. col[%d].j = %d, col[%d].i = %d,
> v[%d]=4\n",procno,j,ncols-1,j,ncols-1,i-1,ncols);
> }
> if (j>0){
> col[ncols].j = j-1; col[ncols].i = i; v[ncols++] = -1.;
> printf("procno: %d, j:%d. col[%d].j = %d, col[%d].i = %d,
> v[%d]=4\n",procno,j,ncols-1,j-1,ncols-1,i,ncols);
> }
> if (j<info.my-1) {
> col[ncols].j = j+1; col[ncols].i = i; v[ncols++] = -1.;
> printf("procno: %d, j:%d. col[%d].j = %d, col[%d].i = %d,
> v[%d]=-1\n",procno,j,ncols-1,j+1,ncols-1,i+1,ncols);
> }
>
> /* boundary: bottom row */
> if (i<info.mx-1) {
> col[ncols].j = j; col[ncols].i = i+1; v[ncols++] = -1.;
> printf("procno: %d, j:%d. col[%d].j = %d, col[%d].i = %d,
> v[%d]=-1\n",procno,j,ncols-1,j,ncols-1,i+1,ncols);
> }
> ierr =
> MatSetValuesStencil(A,1,&row,ncols,col,v,INSERT_VALUES);CHKERRQ(ierr);
> }
> }
>
--
What most experimenters take for granted before they begin their
experiments is infinitely more interesting than any results to which their
experiments lead.
-- Norbert Wiener
https://www.cse.buffalo.edu/~knepley/ <http://www.cse.buffalo.edu/~knepley/>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20220803/0ea4081b/attachment-0001.html>
More information about the petsc-users
mailing list