<html>
<head>
<meta http-equiv="Content-Type" content="text/html; charset=gb2312">
<style type="text/css" style="display:none;"> P {margin-top:0;margin-bottom:0;} </style>
</head>
<body dir="ltr">
<div style="font-family: Calibri, Helvetica, sans-serif; font-size: 12pt; color: rgb(0, 0, 0); background-color: rgb(255, 255, 255);">
Hong,</div>
<div style="font-family: Calibri, Helvetica, sans-serif; font-size: 12pt; color: rgb(0, 0, 0); background-color: rgb(255, 255, 255);">
I added the function to the end, but the results didn't change.</div>
<div style="font-family: Calibri, Helvetica, sans-serif; font-size: 12pt; color: rgb(0, 0, 0); background-color: rgb(255, 255, 255);">
Then I tested another thing, assigning global by Get(Restore)Array instead of DMDANaturalToGlobal. And the results switched. VecView outputs 0 and PetscSynchronizedPrintf outputs correct values. Strange.</div>
<div style="font-family: Calibri, Helvetica, sans-serif; font-size: 12pt; color: rgb(0, 0, 0); background-color: rgb(255, 255, 255);">
<br>
</div>
<div style="font-family: Calibri, Helvetica, sans-serif; font-size: 12pt; color: rgb(0, 0, 0); background-color: rgb(255, 255, 255);">
<pre style="background-color:#2b2b2b;color:#a9b7c6;font-family:'Source Code Pro Semibold';font-size:10.5pt;">    ierr = DMDACreate1d(PETSC_COMM_WORLD<span style="color:#cc7832;">,</span><span style="color:#9876aa;font-style:italic;">DM_BOUNDARY_NONE</span><span style="color:#cc7832;">,</span><span style="color:#6897bb;">5</span><span style="color:#cc7832;">,</span><span style="color:#6897bb;">1</span><span style="color:#cc7832;">,</span><span style="color:#6897bb;">1</span><span style="color:#cc7832;">,</span><span style="color:#908b25;">NULL</span><span style="color:#cc7832;">,</span>&da)<span style="color:#cc7832;">;</span><span style="color:#908b25;">CHKERRQ</span>(ierr)<span style="color:#cc7832;">;<br></span><span style="color:#cc7832;">    </span>ierr = DMSetUp(da)<span style="color:#cc7832;">;</span><span style="color:#908b25;">CHKERRQ</span>(ierr)<span style="color:#cc7832;">;<br></span><span style="color:#cc7832;">    </span>ierr = DMCreateGlobalVector(da<span style="color:#cc7832;">,</span>&global)<span style="color:#cc7832;">;</span><span style="color:#908b25;">CHKERRQ</span>(ierr)<span style="color:#cc7832;">;<br></span><span style="color:#808080;"><br></span><span style="color:#808080;">    </span>ierr = DMDAGetCorners(da<span style="color:#cc7832;">,</span>&x<span style="color:#cc7832;">,</span><span style="color:#908b25;">NULL</span><span style="color:#cc7832;">,</span><span style="color:#908b25;">NULL</span><span style="color:#cc7832;">,</span>&xm<span style="color:#cc7832;">,</span><span style="color:#908b25;">NULL</span><span style="color:#cc7832;">,</span><span style="color:#908b25;">NULL</span>)<span style="color:#cc7832;">;</span><span style="color:#908b25;">CHKERRQ</span>(ierr)<span style="color:#cc7832;">;<br></span><span style="color:#cc7832;">    </span>ierr = DMDAVecGetArray(da<span style="color:#cc7832;">, </span>global<span style="color:#cc7832;">, </span>&val)<span style="color:#cc7832;">;</span><span style="color:#908b25;">CHKERRQ</span>(ierr)<span style="color:#cc7832;">;<br></span><span style="color:#cc7832;">    </span>ierr = PetscSynchronizedPrintf(<span style="color:#908b25;">PETSC_COMM_SELF</span><span style="color:#cc7832;">, </span><span style="color:#6a8759;">"Rank=%d</span><span style="color:#cc7832;">\n</span><span style="color:#6a8759;">"</span><span style="color:#cc7832;">, </span>rank)<span style="color:#cc7832;">;</span><span style="color:#908b25;">CHKERRQ</span>(ierr)<span style="color:#cc7832;">;<br></span><span style="color:#cc7832;">    </span><span style="color:#cc7832;font-weight:bold;">for </span>(i = x<span style="color:#cc7832;">; </span>i < x + xm<span style="color:#cc7832;">; </span>++i) {<br>        val[i] = i<span style="color:#cc7832;">;<br></span><span style="color:#cc7832;">    </span>}<br>    ierr = DMDAVecRestoreArray(da<span style="color:#cc7832;">, </span>global<span style="color:#cc7832;">, </span>&val)<span style="color:#cc7832;">;</span><span style="color:#908b25;">CHKERRQ</span>(ierr)<span style="color:#cc7832;">;<br></span><span style="color:#cc7832;"><br></span><span style="color:#cc7832;">    </span>VecView(global<span style="color:#cc7832;">,</span><span style="color:#908b25;">PETSC_VIEWER_STDOUT_WORLD</span>)<span style="color:#cc7832;">;<br></span><span style="color:#cc7832;">    <br></span><span style="color:#cc7832;">    </span>ierr = DMDAGetCorners(da<span style="color:#cc7832;">,</span>&x<span style="color:#cc7832;">,</span><span style="color:#908b25;">NULL</span><span style="color:#cc7832;">,</span><span style="color:#908b25;">NULL</span><span style="color:#cc7832;">,</span>&xm<span style="color:#cc7832;">,</span><span style="color:#908b25;">NULL</span><span style="color:#cc7832;">,</span><span style="color:#908b25;">NULL</span>)<span style="color:#cc7832;">;</span><span style="color:#908b25;">CHKERRQ</span>(ierr)<span style="color:#cc7832;">;<br></span><span style="color:#cc7832;">    </span>ierr = DMDAVecGetArray(da<span style="color:#cc7832;">, </span>global<span style="color:#cc7832;">, </span>&val)<span style="color:#cc7832;">;</span><span style="color:#908b25;">CHKERRQ</span>(ierr)<span style="color:#cc7832;">;<br></span><span style="color:#cc7832;">    </span>ierr = PetscSynchronizedPrintf(<span style="color:#908b25;">PETSC_COMM_SELF</span><span style="color:#cc7832;">, </span><span style="color:#6a8759;">"Rank=%d</span><span style="color:#cc7832;">\n</span><span style="color:#6a8759;">"</span><span style="color:#cc7832;">, </span>rank)<span style="color:#cc7832;">;</span><span style="color:#908b25;">CHKERRQ</span>(ierr)<span style="color:#cc7832;">;<br></span><span style="color:#cc7832;">    </span><span style="color:#cc7832;font-weight:bold;">for </span>(i = x<span style="color:#cc7832;">; </span>i < x + xm<span style="color:#cc7832;">; </span>++i) {<br>        ierr = PetscSynchronizedPrintf(<span style="color:#908b25;">PETSC_COMM_SELF</span><span style="color:#cc7832;">, </span><span style="color:#6a8759;">"%4.f "</span><span style="color:#cc7832;">, </span>val[i])<span style="color:#cc7832;">;</span><span style="color:#908b25;">CHKERRQ</span>(ierr)<span style="color:#cc7832;">;<br></span><span style="color:#cc7832;">    </span>}<br>    ierr = DMDAVecRestoreArray(da<span style="color:#cc7832;">, </span>global<span style="color:#cc7832;">, </span>&val)<span style="color:#cc7832;">;</span><span style="color:#908b25;">CHKERRQ</span>(ierr)<span style="color:#cc7832;">;<br></span><span style="color:#cc7832;">    </span>PetscSynchronizedFlush(<span style="color:#908b25;">PETSC_COMM_SELF</span><span style="color:#cc7832;">, </span>PETSC_STDOUT)<span style="color:#cc7832;">;</span></pre>
<div>$ mpiexec -n 2 ./test0</div>
<div>Rank=0</div>
<div>Rank=1</div>
<div>Vec Object: 2 MPI processes</div>
<div>  type: mpi</div>
<div>Vec Object: Vec_0x7fffd6629a90_0 2 MPI processes</div>
<div>  type: mpi</div>
<div>Rank=1</div>
<div>   3    4 Process [0]</div>
<div>0.</div>
<div>0.</div>
<div>0.</div>
<div>Process [1]</div>
<div>0.</div>
<div>0.</div>
<div>Rank=0</div>
<div>   0    1    2</div>
<br>
</div>
<div style="font-family: Calibri, Helvetica, sans-serif; font-size: 12pt; color: rgb(0, 0, 0);">
<br>
</div>
<hr style="display:inline-block;width:98%" tabindex="-1">
<div id="divRplyFwdMsg" dir="ltr">
<div> </div>
</div>
<div>
<div dir="ltr">
<div class="x_gmail_extra">
<div class="x_gmail_quote">Songtao :</div>
<div class="x_gmail_quote">You may need adding PetscSynchronizedFlush() to the end. See petsc/src/dm/examples/tests/ex43.c</div>
<div class="x_gmail_quote">Hong<br>
<blockquote class="x_gmail_quote" style="margin:0px 0px 0px 0.8ex; border-left:1px solid rgb(204,204,204); padding-left:1ex">
<div dir="ltr">
<div style="font-family:Calibri,Helvetica,sans-serif; font-size:12pt; color:rgb(0,0,0); background-color:rgb(255,255,255)">
Hello,</div>
<div style="font-family:Calibri,Helvetica,sans-serif; font-size:12pt; color:rgb(0,0,0); background-color:rgb(255,255,255)">
    I am just a beginner. Now I get confused on how to correctly use Vec and DMDA.</div>
<div style="font-family:Calibri,Helvetica,sans-serif; font-size:12pt; color:rgb(0,0,0); background-color:rgb(255,255,255)">
    I set 1D DMDA global vector with natural values for a test. Then use VecView and DMDAVecGetArray to view the data, but the results are different. I don't know why.</div>
<div style="font-family:Calibri,Helvetica,sans-serif; font-size:12pt; color:rgb(0,0,0); background-color:rgb(255,255,255)">
<br>
</div>
<div style="font-family:Calibri,Helvetica,sans-serif; font-size:12pt; color:rgb(0,0,0); background-color:rgb(255,255,255)">
this is the code:</div>
<div style="font-family:Calibri,Helvetica,sans-serif; font-size:12pt; color:rgb(0,0,0); background-color:rgb(255,255,255)">
<pre style="background-color:rgb(43,43,43); color:rgb(169,183,198); font-family:"Source Code Pro Semibold"; font-size:10.5pt"><span style="color:rgb(185,188,209)">PetscMPIInt     </span>rank<span style="color:rgb(204,120,50)">;<br></span><span style="color:rgb(185,188,209)">PetscErrorCode  </span>ierr<span style="color:rgb(204,120,50)">;<br></span><span style="color:rgb(185,188,209)">Vec             </span>global<span style="color:rgb(204,120,50)">,</span>local<span style="color:rgb(204,120,50)">,</span>natural<span style="color:rgb(204,120,50)">;<br></span><span style="color:rgb(185,188,209)">DM              </span>da<span style="color:rgb(204,120,50)">;<br></span><span style="color:rgb(185,188,209)">PetscReal       </span>*xnatural<span style="color:rgb(204,120,50)">,</span>*val<span style="color:rgb(204,120,50)">;<br></span><span style="color:rgb(185,188,209)">PetscInt        </span>i<span style="color:rgb(204,120,50)">,</span>start<span style="color:rgb(204,120,50)">,</span>end<span style="color:rgb(204,120,50)">,</span>x<span style="color:rgb(204,120,50)">,</span>xm<span style="color:rgb(204,120,50)">;<br></span><span style="color:rgb(204,120,50)"><br></span><span style="color:rgb(204,120,50)"><br></span>ierr = PetscInitialize(&argc<span style="color:rgb(204,120,50)">,</span>&argv<span style="color:rgb(204,120,50)">,</span>(<span style="color:rgb(204,120,50); font-weight:bold">c<wbr>har</span>*)<span style="color:rgb(104,151,187)">0</span><span style="color:rgb(204,120,50)">,</span>help)<span style="color:rgb(204,120,50)">;</span><span style="color:rgb(204,120,50); font-weight:bold">if </span>(ierr) <span style="color:rgb(204,120,50); font-weight:bold">return </span>ierr<span style="color:rgb(204,120,50)">;<br></span>ierr = MPI_Comm_rank(PETSC_COMM_WORLD<span style="color:rgb(204,120,50)"><wbr>,</span>&rank)<span style="color:rgb(204,120,50)">;</span><span style="color:rgb(144,139,37)">CHKERRQ</span>(ierr)<span style="color:rgb(204,120,50)">;<br></span><span style="color:rgb(204,120,50)"><br></span>ierr = DMDACreate1d(PETSC_COMM_WORLD<span style="color:rgb(204,120,50)">,</span><span style="color:rgb(152,118,170); font-style:italic"><wbr>DM_BOUNDARY_NONE</span><span style="color:rgb(204,120,50)">,</span><span style="color:rgb(104,151,187)">5</span><span style="color:rgb(204,120,50)">,</span><span style="color:rgb(104,151,187)">1</span><span style="color:rgb(204,120,50)">,</span><span style="color:rgb(104,151,187)">1</span><span style="color:rgb(204,120,50)">,</span><span style="color:rgb(144,139,37)">NULL</span><span style="color:rgb(204,120,50)">,</span>&<wbr>da)<span style="color:rgb(204,120,50)">;</span><span style="color:rgb(144,139,37)">CHKERRQ</span>(ierr)<span style="color:rgb(204,120,50)">;<br></span>ierr = DMSetUp(da)<span style="color:rgb(204,120,50)">;</span><span style="color:rgb(144,139,37)">CHKERRQ</span>(ierr)<span style="color:rgb(204,120,50)">;<br></span>ierr = DMCreateGlobalVector(da<span style="color:rgb(204,120,50)">,</span>&<wbr>global)<span style="color:rgb(204,120,50)">;</span><span style="color:rgb(144,139,37)">CHKERRQ</span>(ierr)<span style="color:rgb(204,120,50)">;<br></span><span style="color:rgb(204,120,50)"><br></span>ierr = DMDACreateNaturalVector(da<span style="color:rgb(204,120,50)">,</span>&<wbr>natural)<span style="color:rgb(204,120,50)">;</span><span style="color:rgb(144,139,37)">CHKERRQ</span>(ierr)<span style="color:rgb(204,120,50)">;<br></span>ierr = VecGetOwnershipRange(natural<span style="color:rgb(204,120,50)">,</span>&<wbr>start<span style="color:rgb(204,120,50)">,</span>&end)<span style="color:rgb(204,120,50)">;</span><span style="color:rgb(144,139,37)">CHKERRQ</span>(ierr)<span style="color:rgb(204,120,50)">;<br></span>ierr = VecGetArray(natural<span style="color:rgb(204,120,50)">,</span>&xnatural)<span style="color:rgb(204,120,50)"><wbr>;</span><span style="color:rgb(144,139,37)">CHKERRQ</span>(ierr)<span style="color:rgb(204,120,50)">;<br></span><span style="color:rgb(204,120,50); font-weight:bold">for </span>(i=start<span style="color:rgb(204,120,50)">; </span>i<end<span style="color:rgb(204,120,50)">; </span>i++) {<br>    xnatural[i-start] = i<span style="color:rgb(204,120,50)">;<br></span>}<br>ierr = VecRestoreArray(natural<span style="color:rgb(204,120,50)">,</span>&<wbr>xnatural)<span style="color:rgb(204,120,50)">;</span><span style="color:rgb(144,139,37)">CHKERRQ</span>(ierr)<span style="color:rgb(204,120,50)">;<br></span>ierr = DMDANaturalToGlobalBegin(da<span style="color:rgb(204,120,50)">,</span>na<wbr>tural<span style="color:rgb(204,120,50)">,</span><span style="color:rgb(152,118,170); font-style:italic">INSERT_VALUES</span><span style="color:rgb(204,120,50)">,</span>global)<span style="color:rgb(204,120,50)">;</span><span style="color:rgb(144,139,37)">CH<wbr>KERRQ</span>(ierr)<span style="color:rgb(204,120,50)">;<br></span>ierr = DMDANaturalToGlobalEnd(da<span style="color:rgb(204,120,50)">,</span>natu<wbr>ral<span style="color:rgb(204,120,50)">,</span><span style="color:rgb(152,118,170); font-style:italic">INSERT_VALUES</span><span style="color:rgb(204,120,50)">,</span>global)<span style="color:rgb(204,120,50)">;</span><span style="color:rgb(144,139,37)">CHKE<wbr>RRQ</span>(ierr)<span style="color:rgb(204,120,50)">;<br></span>ierr = VecDestroy(&natural)<span style="color:rgb(204,120,50)">;</span><span style="color:rgb(144,139,37)">CHKERRQ</span>(<wbr>ierr)<span style="color:rgb(204,120,50)">;<br></span><span style="color:rgb(204,120,50)"><br></span>VecView(global<span style="color:rgb(204,120,50)">,</span><span style="color:rgb(144,139,37)">PETSC_VIEWER_<wbr>STDOUT_WORLD</span>)<span style="color:rgb(204,120,50)">;<br></span><span style="color:rgb(204,120,50)"><br></span>ierr = DMDAGetCorners(da<span style="color:rgb(204,120,50)">,</span>&x<span style="color:rgb(204,120,50)">,</span><span style="color:rgb(144,139,37)">NULL</span><span style="color:rgb(204,120,50)">,</span><span style="color:rgb(144,139,37)">NULL</span><span style="color:rgb(204,120,50)"><wbr>,</span>&xm<span style="color:rgb(204,120,50)">,</span><span style="color:rgb(144,139,37)">NULL</span><span style="color:rgb(204,120,50)">,</span><span style="color:rgb(144,139,37)">NULL</span>)<span style="color:rgb(204,120,50)">;</span><span style="color:rgb(144,139,37)">CHKERRQ</span>(ierr)<span style="color:rgb(204,120,50)">;<br></span>ierr = DMDAVecGetArray(da<span style="color:rgb(204,120,50)">, </span>global<span style="color:rgb(204,120,50)">, </span>&val)<span style="color:rgb(204,120,50)">;</span><span style="color:rgb(144,139,37)">CHKERRQ</span>(ierr)<span style="color:rgb(204,120,50)">;<br></span>ierr = PetscSynchronizedPrintf(<span style="color:rgb(144,139,37)">PETSC_<wbr>COMM_SELF</span><span style="color:rgb(204,120,50)">, </span><span style="color:rgb(106,135,89)">"Rank=%d</span><span style="color:rgb(204,120,50)">\n</span><span style="color:rgb(106,135,89)">"</span><span style="color:rgb(204,120,50)">, </span>rank)<span style="color:rgb(204,120,50)">;</span><span style="color:rgb(144,139,37)">CHKERRQ</span>(ierr)<span style="color:rgb(204,120,50)">;<br></span><span style="color:rgb(204,120,50); font-weight:bold">for </span>(i = x<span style="color:rgb(204,120,50)">; </span>i < x + xm<span style="color:rgb(204,120,50)">; </span>++i) {<br>    ierr = PetscSynchronizedPrintf(<span style="color:rgb(144,139,37)">PETSC_<wbr>COMM_SELF</span><span style="color:rgb(204,120,50)">, </span><span style="color:rgb(106,135,89)">"%4.f "</span><span style="color:rgb(204,120,50)">, </span>val[i])<span style="color:rgb(204,120,50)">;</span><span style="color:rgb(144,139,37)">CHKERRQ</span>(ierr)<span style="color:rgb(204,120,50)">;<br></span>}<br>ierr = DMDAVecRestoreArray(da<span style="color:rgb(204,120,50)">, </span>global<span style="color:rgb(204,120,50)">, </span>&val)<span style="color:rgb(204,120,50)">;</span><span style="color:rgb(144,139,37)">CHKERRQ</span>(ierr)<span style="color:rgb(204,120,50)">;<br></span></pre>
and the results are:</div>
<div style="font-family:Calibri,Helvetica,sans-serif; font-size:12pt; color:rgb(0,0,0); background-color:rgb(255,255,255)">
<br>
</div>
<div style="font-family:Calibri,Helvetica,sans-serif; font-size:12pt; color:rgb(0,0,0); background-color:rgb(255,255,255)">
<div>$ ./test0</div>
<div>Vec Object: 1 MPI processes</div>
<div>  type: seq</div>
<div>Vec Object: Vec_0x7ffff3cbc500_0 1 MPI processes</div>
<div>  type: mpi</div>
<div>Process [0]</div>
<div>0.</div>
<div>1.</div>
<div>2.</div>
<div>3.</div>
<div>4.</div>
<div>Rank=0</div>
<div>   0    1    2    3    4 </div>
<div><br>
</div>
<div>$ mpiexec -n 2 ./test0</div>
<div>Vec Object: 2 MPI processes</div>
<div>  type: mpi</div>
<div>Vec Object: Vec_0x7fffcf948a90_0 2 MPI processes</div>
<div>  type: mpi</div>
<div>Process [0]</div>
<div>0.</div>
<div>1.</div>
<div>2.</div>
<div>Process [1]</div>
<div>3.</div>
<div>4.</div>
<div>Rank=0</div>
<div>   0    0    0 Rank=1</div>
<div>   0    0</div>
<br>
</div>
<div style="font-family:Calibri,Helvetica,sans-serif; font-size:12pt; color:rgb(0,0,0); background-color:rgb(255,255,255)">
<br>
</div>
<div style="font-family:Calibri,Helvetica,sans-serif; font-size:12pt; color:rgb(0,0,0); background-color:rgb(255,255,255)">
the seq version is correct, but the mpi version is wrong. Every element in val is 0.</div>
<div style="font-family:Calibri,Helvetica,sans-serif; font-size:12pt; color:rgb(0,0,0); background-color:rgb(255,255,255)">
<br>
</div>
</div>
</blockquote>
</div>
<br>
</div>
</div>
</div>
</body>
</html>