On Tue, Apr 3, 2012 at 7:40 AM, Barry Smith <span dir="ltr">&lt;<a href="mailto:bsmith@mcs.anl.gov">bsmith@mcs.anl.gov</a>&gt;</span> wrote:<br><div class="gmail_quote"><blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex">
<br>
   I think you are mis-understanding the global ordering when using DA vectors in parallel. You might consider using VecSetValuesLocal() or VecSetValuesStencil() instead of VecSetValues(). Read up on the users manual and examples about PETSc global ordering.</blockquote>
<div><br></div><div>More explicitly, you do set the 5th value of the vector each time, but the layout is different when created with the DA</div><div>since it is splitting the grid somehow. When you call VecView() it reorders the Vec to match the grid before output.</div>
<div><br></div><div>    Matt</div><div> </div><blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex"><span class="HOEnZb"><font color="#888888"><br>
   Barry<br>
</font></span><div class="HOEnZb"><div class="h5"><br>
On Apr 3, 2012, at 3:46 AM, khalid ashraf wrote:<br>
<br>
&gt; Hi Matt, I ran the program on two processors. I am trying to enter value into a<br>
&gt; vector associated with DA. I am using VecSetValues since I may need to assign nonlocal<br>
&gt; values to the vector. However, I find that on 2 procs, a general parallel vector works but<br>
&gt; the DAVector indexing is not working. Do you see any problem with the DA code ?<br>
&gt;<br>
&gt; Thanks.<br>
&gt;<br>
&gt; Khalid<br>
&gt;<br>
&gt; THIS WORKS: 30 is written at index 5<br>
&gt; -------------------------------<br>
&gt; ppy=5; ppx=30;<br>
&gt; VecCreate(PETSC_COMM_WORLD,&amp;testVec);<br>
&gt; VecSetFromOptions(testVec);<br>
&gt; VecSetSizes(testVec,PETSC_DECIDE,384);<br>
&gt;<br>
&gt; VecSetValues(testVec,1,&amp;ppy,&amp;ppx,INSERT_VALUES);<br>
&gt; VecAssemblyBegin(testVec);<br>
&gt; VecAssemblyEnd(testVec);<br>
&gt; =====================================<br>
&gt;<br>
&gt; THIS DOES NOT WORK: puts 30 at index 9 instead of 5 when run on two procs.<br>
&gt; ---------------------------------<br>
&gt; ppy=5; ppx=30;<br>
&gt; ierr = DACreate3d(PETSC_COMM_WORLD,DA_YPERIODIC,DA_STENCIL_BOX,8,8,6,<br>
&gt;                     PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,1,1,PETSC_NULL,PETSC_NULL,PETSC_NULL,<br>
&gt;                     &amp;appctx.da);CHKERRQ(ierr);<br>
&gt;  ierr = DACreateGlobalVector(da,&amp;testVec);CHKERRQ(ierr)<br>
&gt; VecSetValues(testVec,1,&amp;ppy,&amp;ppx,INSERT_VALUES);<br>
&gt; VecAssemblyBegin(testVec);<br>
&gt; VecAssemblyEnd(testVec);<br>
&gt;<br>
</div></div><div class="HOEnZb"><div class="h5">&gt; From: Matthew Knepley &lt;<a href="mailto:knepley@gmail.com">knepley@gmail.com</a>&gt;<br>
&gt; To: khalid ashraf &lt;<a href="mailto:khalid_eee@yahoo.com">khalid_eee@yahoo.com</a>&gt;<br>
&gt; Cc: PETSc users list &lt;<a href="mailto:petsc-users@mcs.anl.gov">petsc-users@mcs.anl.gov</a>&gt;<br>
&gt; Sent: Monday, April 2, 2012 3:08 PM<br>
&gt; Subject: Re: [petsc-users] transfer vector data diagonally on DA<br>
&gt;<br>
&gt; On Mon, Apr 2, 2012 at 4:53 PM, khalid ashraf &lt;<a href="mailto:khalid_eee@yahoo.com">khalid_eee@yahoo.com</a>&gt; wrote:<br>
&gt; Thanks Matt. It works now.<br>
&gt;<br>
&gt; I have another problem. I am writing this program.<br>
&gt;<br>
&gt; 1) Always test on 2 procs once<br>
&gt;<br>
&gt; 2) I have no idea what you actually want to do. However, it really simple to just print out<br>
&gt;     the indices you are using in row[] and the values in bodyFx[] and check that VecSetValues()<br>
&gt;     is working as you expect.<br>
&gt;<br>
&gt;    Matt<br>
&gt;<br>
&gt;<br>
&gt;  for (k=zs; k&lt;zs+zm; k++) {<br>
&gt;   for (j=ys; j&lt;ys+ym; j++) {<br>
&gt;   for (i=xs; i&lt;xs+xm; i++) {<br>
&gt;  if ( i!=(mx-1) || j!=my-1 || k!=mz-1)<br>
&gt; {<br>
&gt; bodyFx[0]=1;<br>
&gt; bodyFx[1]=-0.3;<br>
&gt; bodyFx[2]=-0.3;<br>
&gt; bodyFx[3]=-0.3;<br>
&gt; row[0]=k*mx*my+j*mx+i;<br>
&gt; row[1]=k*mx*my+j*mx+i+1;<br>
&gt; row[2]=k*mx*my+(j+1)*mx+i;<br>
&gt; row[3]=(k+1)*mx*my+j*mx+i;<br>
&gt; VecSetValues(fx_test,4,row,bodyFx,ADD_VALUES);<br>
&gt; }<br>
&gt; VecAssemblyBegin(fx_test);<br>
&gt; VecAssemblyEnd(fx_test);<br>
&gt;<br>
&gt; Output: Print fx_test<br>
&gt;<br>
&gt;  Here fx_test is a global vector. This program gives correct output on 1 proc. but on 4 processors,<br>
&gt; the program runs but the outputs are wrong.<br>
&gt;<br>
&gt; I also tried<br>
&gt;  for (k=zs; k&lt;zs+zm; k++) {<br>
&gt;   for (j=ys; j&lt;ys+ym; j++) {<br>
&gt;   for (i=xs; i&lt;xs+xm; i++) {<br>
&gt;  if ( i!=(mx-1) || j!=my-1 || k!=mz-1)<br>
&gt; {<br>
&gt; bodyFx[0]=1;<br>
&gt; bodyFx[1]=-0.3;<br>
&gt; bodyFx[2]=-0.3;<br>
&gt; bodyFx[3]=-0.3;<br>
&gt;  fx1_localptr[k][j][i+1]+=bodyFx[4];<br>
&gt;  fx1_localptr[k][j+1][i]+=bodyFx[3];<br>
&gt;  fx1_localptr[k+1][j][i]+=bodyFx[1];<br>
&gt;  fx1_localptr[k][j][i]+=bodyFx[0];<br>
&gt; }<br>
&gt; Here, fx1_localptr is an array from a local vector that is communicated from the global vector fx_test.<br>
&gt; In this case as well, I get similar errors.<br>
&gt;<br>
&gt; Thanks.<br>
&gt; Khalid<br>
&gt;<br>
&gt; Some final numbers of the vector with 8*8*6 grid, 4 proc,<br>
&gt; -2.9999999999999999e-01<br>
&gt; -2.9999999999999999e-01<br>
&gt; -2.9999999999999999e-01<br>
&gt; -2.9999999999999999e-01<br>
&gt; 1.0000000000000009e-01<br>
&gt; 1.0000000000000009e-01<br>
&gt; 1.0000000000000009e-01<br>
&gt; -2.9999999999999999e-01<br>
&gt; -2.9999999999999999e-01<br>
&gt; -2.9999999999999999e-01<br>
&gt; -2.9999999999999999e-01<br>
&gt; 0.0000000000000000e+00<br>
&gt; 4.0000000000000002e-01<br>
&gt; 1.0000000000000009e-01<br>
&gt; 1.0000000000000009e-01<br>
&gt; 1.0000000000000009e-01<br>
&gt; 0.0000000000000000e+00<br>
&gt; 0.0000000000000000e+00<br>
&gt; 0.0000000000000000e+00<br>
&gt; 0.0000000000000000e+00<br>
&gt; 1.0000000000000009e-01<br>
&gt; 1.0000000000000009e-01<br>
&gt; 1.0000000000000009e-01<br>
&gt; -2.9999999999999999e-01<br>
&gt; 0.0000000000000000e+00<br>
&gt; 0.0000000000000000e+00<br>
&gt; 0.0000000000000000e+00<br>
&gt;<br>
&gt;<br>
&gt; Some final numbers of the vector with 8*8*6 grid, 1 proc,<br>
&gt; 0.0000000000000000e+00<br>
&gt; -2.9999999999999999e-01<br>
&gt; -2.9999999999999999e-01<br>
&gt; -2.9999999999999999e-01<br>
&gt; -2.9999999999999999e-01<br>
&gt; -2.9999999999999999e-01<br>
&gt; -2.9999999999999999e-01<br>
&gt; -2.9999999999999999e-01<br>
&gt; 0.0000000000000000e+00<br>
&gt; -2.9999999999999999e-01<br>
&gt; -2.9999999999999999e-01<br>
&gt; -2.9999999999999999e-01<br>
&gt; -2.9999999999999999e-01<br>
&gt; -2.9999999999999999e-01<br>
&gt; -2.9999999999999999e-01<br>
&gt; -2.9999999999999999e-01<br>
&gt; 0.0000000000000000e+00<br>
&gt; -2.9999999999999999e-01<br>
&gt; -2.9999999999999999e-01<br>
&gt; -2.9999999999999999e-01<br>
&gt; -2.9999999999999999e-01<br>
&gt; -2.9999999999999999e-01<br>
&gt; -2.9999999999999999e-01<br>
&gt; -2.9999999999999999e-01<br>
&gt; 0.0000000000000000e+00<br>
&gt; 0.0000000000000000e+00<br>
&gt; 0.0000000000000000e+00<br>
&gt; 0.0000000000000000e+00<br>
&gt; 0.0000000000000000e+00<br>
&gt; 0.0000000000000000e+00<br>
&gt; 0.0000000000000000e+00<br>
&gt; 0.0000000000000000e+00<br>
&gt; 0.0000000000000000e+00<br>
&gt;<br>
&gt;<br>
&gt;<br>
</div></div><div class="HOEnZb"><div class="h5">&gt; From: Matthew Knepley &lt;<a href="mailto:knepley@gmail.com">knepley@gmail.com</a>&gt;<br>
&gt; To: khalid ashraf &lt;<a href="mailto:khalid_eee@yahoo.com">khalid_eee@yahoo.com</a>&gt;; PETSc users list &lt;<a href="mailto:petsc-users@mcs.anl.gov">petsc-users@mcs.anl.gov</a>&gt;<br>
&gt; Sent: Monday, April 2, 2012 5:41 AM<br>
&gt; Subject: Re: [petsc-users] transfer vector data diagonally on DA<br>
&gt;<br>
&gt; On Mon, Apr 2, 2012 at 3:36 AM, khalid ashraf &lt;<a href="mailto:khalid_eee@yahoo.com">khalid_eee@yahoo.com</a>&gt; wrote:<br>
&gt; Hi Jed,<br>
&gt;<br>
&gt; I am using petsc/3.1 and the include file &quot;petscdmda.h&quot; is not working. I am using the &quot;petscda.h&quot;<br>
&gt; and DACreate3D.<br>
&gt; unew_localptr[][][] is from a global vector using VecGetArray<br>
&gt; u_localptr[][][] is from a local vector communicated from a global vector.<br>
&gt; I tried<br>
&gt; unew_localptr[k][j][i]=u_localptr[k-1][j-1][i-1]<br>
&gt; but it gives Segmentation Violation error.<br>
&gt;<br>
&gt; I would guess (because you provide almost no information about what you are doing), that<br>
&gt; this is a domain on the edge. If you truly mean to have a periodic domain, you must set that in the<br>
&gt; creation call. Then the local vector will also ghost regions outside the domain boundary.<br>
&gt;<br>
&gt;    Matt<br>
&gt;<br>
&gt; I also tried unew_localptr[k][j][i]=u_localptr[k][j+1][i+1] which works but<br>
&gt; unew_localptr[k][j][i]=u_localptr[k+1][j+1][i+1] or unew_localptr[k][j][i]=u_localptr[k][j-1][i-1]<br>
&gt; does not work. I am running the program on4 processors.<br>
&gt;<br>
&gt; Thanks.<br>
&gt;<br>
&gt;<br>
&gt;   ierr = DACreate3d(PETSC_COMM_WORLD,DA_YPERIODIC,DA_STENCIL_BOX,appctx.l,appctx.m,appctx.n,<br>
&gt;                     PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,1,1,PETSC_NULL,PETSC_NULL,PETSC_NULL,<br>
&gt;                     &amp;appctx.da);CHKERRQ(ierr);<br>
&gt;  for (k=zs; k&lt;zs+zm; k++) {<br>
&gt;   for (j=ys; j&lt;ys+ym; j++) {<br>
&gt;   for (i=xs; i&lt;xs+xm; i++) {<br>
&gt;<br>
&gt;   if(i!=0 || j!=0 || k!=0|| i!=mx-1 || j!=my-1 || k!=mz-1)<br>
&gt;<br>
&gt;   unew_localptr[k][j][i]=u_localptr[k+1][j+1][i+1];<br>
&gt; }}}<br>
</div></div><div class="HOEnZb"><div class="h5">&gt; From: Jed Brown &lt;<a href="mailto:jedbrown@mcs.anl.gov">jedbrown@mcs.anl.gov</a>&gt;<br>
&gt; To: khalid ashraf &lt;<a href="mailto:khalid_eee@yahoo.com">khalid_eee@yahoo.com</a>&gt;; PETSc users list &lt;<a href="mailto:petsc-users@mcs.anl.gov">petsc-users@mcs.anl.gov</a>&gt;<br>
&gt; Sent: Sunday, April 1, 2012 10:07 PM<br>
&gt; Subject: Re: [petsc-users] transfer vector data diagonally on DA<br>
&gt;<br>
&gt; On Sun, Apr 1, 2012 at 22:01, khalid ashraf &lt;<a href="mailto:khalid_eee@yahoo.com">khalid_eee@yahoo.com</a>&gt; wrote:<br>
&gt; I want to transfer vector data diagonally in the DA grid like<br>
&gt;  for (k=zs; k&lt;zs+zm; k++) {<br>
&gt;   for (j=ys; j&lt;ys+ym; j++) {<br>
&gt;   for (i=xs; i&lt;xs+xm; i++) {<br>
&gt; if(i!=mx-1 || j!=my-1 || k!=mz-1){<br>
&gt; u_new[k+1][j+1][i+1]=u[k][j][i];}<br>
&gt; }}}<br>
&gt;<br>
&gt; Could you please suggest the best way to do it minimizing interprocessor assignments.<br>
&gt;<br>
&gt; Both are on the same DMDA?<br>
&gt;<br>
&gt; Communicate U to Ulocal (DMGlobalToLocalBegin/End) using a BOX stencil with width at least 1, get the global array u_new[][][] from UGlobalNew and the local arrays u[][][] from Ulocal, then assign u_new[k][j][i] = u[k-1][j-1][i-1].<br>

&gt;<br>
&gt;<br>
&gt;<br>
&gt;<br>
&gt;<br>
&gt; --<br>
&gt; What most experimenters take for granted before they begin their experiments is infinitely more interesting than any results to which their experiments lead.<br>
&gt; -- Norbert Wiener<br>
&gt;<br>
&gt;<br>
&gt;<br>
&gt;<br>
&gt;<br>
&gt; --<br>
&gt; What most experimenters take for granted before they begin their experiments is infinitely more interesting than any results to which their experiments lead.<br>
&gt; -- Norbert Wiener<br>
&gt;<br>
&gt;<br>
<br>
</div></div></blockquote></div><br><br clear="all"><div><br></div>-- <br>What most experimenters take for granted before they begin their experiments is infinitely more interesting than any results to which their experiments lead.<br>
-- Norbert Wiener<br>