<html><head><meta http-equiv="Content-Type" content="text/html; charset=us-ascii"></head><body style="word-wrap: break-word; -webkit-nbsp-mode: space; line-break: after-white-space;" class=""><div class=""><br class=""></div> What does "not working as I would like" mean? It should be retrieving the trajectory at the times 1.0, 2.0, 3.0 ... 40.0 and setting into the vector partial the values of the second component of Uloc (which depending on DMDA having a stencil width of 1 and a w of 1 is the first component of U.<div class=""><br class=""></div><div class=""> You can move the <span style="font-family: Menlo; font-size: 11px;" class="">VecGet/RestoreArray(partial,&partlocal);CHKERRQ(ierr); </span>outside of the loop.</div><div class=""><br class=""></div><div class=""> If you want the first component of U on process 0 you don't need the Uloc or the GlobalToLocalBegin/End. just use <span style="font-family: Menlo; font-size: 11px;" class="">DMDAVecGetArray(appctx.da,U,&Ui);CHKERRQ(ierr);</span></div><div class=""><font face="Menlo" class=""><span style="font-size: 11px;" class=""><br class=""></span></font></div><div class=""><font face="Menlo" class=""><span style="font-size: 11px;" class=""> <br class=""></span></font><div> You only provide 14 locations in partial distributed over the MPI ranks but likely you want 40 on the first rank and none on the other ranks</div><div><br class=""></div><div> You are assigning part local[i] on all ranks, but you said you only want it on rank 0 so here is code that may work</div><div><br class=""></div><div> if rank == 0 {</div><div><blockquote type="cite" class=""><div class="" style="word-wrap: break-word; -webkit-nbsp-mode: space; line-break: after-white-space;"><div class=""><span class="" style="font-size: 11px;"><font face="Menlo" class="">ierr = VecCreateMPI(PETSC_COMM_WORLD,40,PETSC_DETERMINE,&partial);CHKERRQ(ierr); /* 40 local values</font></span></div></div></blockquote><blockquote type="cite" class=""><div class="" style="word-wrap: break-word; -webkit-nbsp-mode: space; line-break: after-white-space;"><div class=""><span class="" style="font-size: 11px;"><font face="Menlo" class=""> ierr = VecGetArray(partial,&partlocal);CHKERRQ(ierr);</font></span></div></div></blockquote>} else {</div><div>ierr = VecCreateMPI(PETSC_COMM_WORLD,0,PETSC_DETERMINE,&partial);CHKERRQ(ierr); /* 0 local values</div><div>}<br class=""><blockquote type="cite" class=""><div class="" style="word-wrap: break-word; -webkit-nbsp-mode: space; line-break: after-white-space;"><div class=""><span class="" style="font-size: 11px;"><font face="Menlo" class=""> for(i=0; i<40; i++) {<br class=""> PetscReal ttime = i+1;<br class=""> ierr = TSTrajectoryGetVecs(appctx.tj,appctx.ts,PETSC_DECIDE,&ttime,U,NULL);CHKERRQ(ierr);<br class=""></font></span></div></div></blockquote> if rank == 0 {<br class=""><blockquote type="cite" class=""><div class="" style="word-wrap: break-word; -webkit-nbsp-mode: space; line-break: after-white-space;"><div class=""><span class="" style="font-size: 11px;"><font face="Menlo" class=""> ierr = DMDAVecGetArray(appctx.da,U,&Ui);CHKERRQ(ierr);<br class=""> partlocal[i] = Ui[0];<br class=""> ierr = DMDAVecRestoreArray(appctx.da,U,&Ui);CHKERRQ(ierr);<br class=""></font></span></div></div></blockquote> }</div><div> }<br class=""><blockquote type="cite" class=""><div class="" style="word-wrap: break-word; -webkit-nbsp-mode: space; line-break: after-white-space;"><div class=""><span class="" style="font-size: 11px;"><font face="Menlo" class=""> if rank == 0 {</font></span></div></div></blockquote><blockquote type="cite" class=""><div class="" style="word-wrap: break-word; -webkit-nbsp-mode: space; line-break: after-white-space;"><div class=""><span class="" style="font-size: 11px;"><font face="Menlo" class=""> ierr = VecRestoreArray(partial,&partlocal);CHKERRQ(ierr);</font></span></div></div></blockquote> }<br class=""><div><br class=""></div> Note that this entire block of code needs to run on all MPI ranks. But the actually selection of the wanted value only occurs on rank 0</div><div><br class=""></div><div> When the loop is done rank == 0 will have a parallel vector whose components are what you want and all the other ranks will have a parallel vector with no components on those ranks. Note that you don't need to make partial be a parallel vector, you can just make it live on rank == 0 because that is the only place you access it. Then the code would be simpler</div><div><br class=""></div><div><div>if rank == 0 {</div><div><blockquote type="cite" class=""><div class="" style="word-wrap: break-word; -webkit-nbsp-mode: space; line-break: after-white-space;"><div class=""><span class="" style="font-size: 11px;"><font face="Menlo" class=""> ierr = VecCreateSeq(PETSC_COMM_WORLD,40,PETSC_&partial);CHKERRQ(ierr);</font></span></div></div></blockquote><blockquote type="cite" class=""><div class="" style="word-wrap: break-word; -webkit-nbsp-mode: space; line-break: after-white-space;"><div class=""><span class="" style="font-size: 11px;"><font face="Menlo" class=""> ierr = VecGetArray(partial,&partlocal);CHKERRQ(ierr);</font></span></div></div></blockquote>}</div><div><span style="font-family: Menlo; font-size: 11px;" class=""> </span></div><div><span style="font-family: Menlo; font-size: 11px;" class="">for(i=0; i<40; i++) {</span><br class=""><blockquote type="cite" class=""><div class="" style="word-wrap: break-word; -webkit-nbsp-mode: space; line-break: after-white-space;"><div class=""><span class="" style="font-size: 11px;"><font face="Menlo" class=""> PetscReal ttime = i+1;<br class=""> ierr = TSTrajectoryGetVecs(appctx.tj,appctx.ts,PETSC_DECIDE,&ttime,U,NULL);CHKERRQ(ierr);<br class=""></font></span></div></div></blockquote> if rank == 0 {<br class=""><blockquote type="cite" class=""><div class="" style="word-wrap: break-word; -webkit-nbsp-mode: space; line-break: after-white-space;"><div class=""><span class="" style="font-size: 11px;"><font face="Menlo" class=""> ierr = DMDAVecGetArray(appctx.da,U,&Ui);CHKERRQ(ierr);<br class=""> partlocal[i] = Ui[0];<br class=""> ierr = DMDAVecRestoreArray(appctx.da,U,&Ui);CHKERRQ(ierr);<br class=""></font></span></div></div></blockquote> }</div><div> }<br class=""><blockquote type="cite" class=""><div class="" style="word-wrap: break-word; -webkit-nbsp-mode: space; line-break: after-white-space;"><div class=""><span class="" style="font-size: 11px;"><font face="Menlo" class=""> if rank == 0 {</font></span></div></div></blockquote><blockquote type="cite" class=""><div class="" style="word-wrap: break-word; -webkit-nbsp-mode: space; line-break: after-white-space;"><div class=""><span class="" style="font-size: 11px;"><font face="Menlo" class=""> ierr = VecRestoreArray(partial,&partlocal);CHKERRQ(ierr);}</font></span></div></div></blockquote><div class=""><div class="" style="word-wrap: break-word; -webkit-nbsp-mode: space; line-break: after-white-space;"><div class=""><span class="" style="font-size: 11px;"><font face="Menlo" class=""><br class=""></font></span></div><div class=""><span class="" style="font-size: 11px;"><font face="Menlo" class=""><br class=""></font></span></div></div></div></div></div><div> Barry</div><div><br class=""><br class=""><blockquote type="cite" class=""><div class="">On May 27, 2021, at 2:42 AM, Francesco Brarda <<a href="mailto:brardafrancesco@gmail.com" class="">brardafrancesco@gmail.com</a>> wrote:</div><br class="Apple-interchange-newline"><div class=""><meta http-equiv="Content-Type" content="text/html; charset=us-ascii" class=""><div style="word-wrap: break-word; -webkit-nbsp-mode: space; line-break: after-white-space;" class="">I created a for cycle where I call TSTrajectoryGetVecs, but only the 0 rank seems to enter in this cycle and I do not know why.<div class="">I thought the following might be a solution, but it is not working as I would like to, i.e. the final vector has the same local parts, a copy of the values obtained with the 0-rank. How should I change this, please?</div><div class=""><div class=""><br class=""></div><span style="font-size: 11px;" class=""><font face="Menlo" class=""><font class="">Vec U, partial, Uloc;<br class="">PetscScalar *Ui, *partlocal;<br class="">PetscInt i;</font><br class="">ierr = VecCreateMPI(PETSC_COMM_WORLD,PETSC_DECIDE,14,&partial);CHKERRQ(ierr);<br class=""> for(i=0; i<40; i++) {<br class=""> PetscReal ttime = i+1;<br class=""> ierr = TSTrajectoryGetVecs(appctx.tj,appctx.ts,PETSC_DECIDE,&ttime,U,NULL);CHKERRQ(ierr);<br class=""> ierr = DMGetLocalVector(appctx.da,&Uloc);CHKERRQ(ierr);<br class=""> ierr = DMGlobalToLocalBegin(appctx.da,U,INSERT_VALUES,Uloc);CHKERRQ(ierr);<br class=""> ierr = DMGlobalToLocalEnd(appctx.da,U,INSERT_VALUES,Uloc);CHKERRQ(ierr);<br class=""> ierr = DMDAVecGetArray(appctx.da,Uloc,&Ui);CHKERRQ(ierr);<br class=""> ierr = VecGetArray(partial,&partlocal);CHKERRQ(ierr);<br class=""> partlocal[i] = Ui[1];<br class=""> ierr = DMDAVecRestoreArray(appctx.da,U,&Ui);CHKERRQ(ierr);<br class=""> ierr = VecRestoreArray(partial,&partlocal);CHKERRQ(ierr);<br class=""> ierr = DMRestoreLocalVector(appctx.da,&Uloc);CHKERRQ(ierr); </font></span></div><div class=""><span class=""><font face="Menlo" style="font-size: 11px;" class="">}</font><br class=""><br class=""><br class=""><blockquote type="cite" class="">Il giorno 27 mag 2021, alle ore 01:15, Barry Smith <<a href="mailto:bsmith@petsc.dev" class="">bsmith@petsc.dev</a>> ha scritto:<br class=""><br class=""><br class=""><br class=""><blockquote type="cite" class="">On May 26, 2021, at 10:39 AM, Francesco Brarda <<a href="mailto:brardafrancesco@gmail.com" class="">brardafrancesco@gmail.com</a>> wrote:<br class=""><br class="">Thank you very much.<br class=""><blockquote type="cite" class=""> Based on the error message it appears that your code is requesting different times on different MPI ranks. Is that what you intend to do?<br class=""></blockquote>Yes. I want to save different times across a vector built with multiple MPI ranks (PETSC_DECIDE for the local length).<br class="">The function is called only by the first proc (rank=0) and not from the others. Is there a way to force also other ranks to call that routine?<br class=""></blockquote><br class=""> Yes, just have all ranks call it and ignore the result on the other ranks. <br class=""><br class=""><blockquote type="cite" class="">Should I build everything into an external function outside the main?<br class=""></blockquote><br class=""> It can be called in main, does not need to be in a different function.<br class=""><br class=""><blockquote type="cite" class=""><br class="">Francesco<br class=""><br class=""><br class=""><blockquote type="cite" class="">Il giorno 26 mag 2021, alle ore 16:20, Barry Smith <<a href="mailto:bsmith@petsc.dev" class="">bsmith@petsc.dev</a>> ha scritto:<br class=""><br class=""><br class=""> <br class=""><br class=""> TSTrajectoryGetVecs() is listed as Collective on TS. This means all ranks must call it with the same times in the same order of operations on all ranks that share the TS. <br class=""><br class=""> You do not need to use VecScatter. Each process must call TSTrajectoryGetVecs with the same time but then you can have only the rank you care about select the entries from the resulting vectors you care about while the other ranks for that time just ignore the vectors since they do not need to values from it. <br class=""><br class=""> Barry<br class=""><br class=""> <br class=""><br class=""><blockquote type="cite" class="">On May 26, 2021, at 5:20 AM, Francesco Brarda <<a href="mailto:brardafrancesco@gmail.com" class="">brardafrancesco@gmail.com</a>> wrote:<br class=""><br class="">Hi!<br class=""><br class="">I solved an ODE system with TS. Now I would like to save one of the trajectories in specific times. To do so, I used TSTrajectoryGetVecs. <br class="">The values of the variable I am interested in is on one processor. I want to collect these values in a parallel vector, but I had the error:<br class=""><br class="">[0]PETSC ERROR: Invalid argument<br class="">[0]PETSC ERROR: Real value must be same on all processes, argument # 2<br class="">[0]PETSC ERROR: See <a href="https://www.mcs.anl.gov/petsc/documentation/faq.html" class="">https://www.mcs.anl.gov/petsc/documentation/faq.html</a> for trouble shooting.<br class="">[0]PETSC ERROR: Petsc Release Version 3.14.4, unknown <br class="">[0]PETSC ERROR: ./petsc_sir on a arch-debug named srvulx13 by fbrarda Wed May 26 12:00:42 2021<br class="">[0]PETSC ERROR: Configure options --with-cc=gcc --with-cxx=g++ --with-fc=gfortran --with-openblas-dir=/opt/packages/openblas/0.2.13-gcc --download-mpich PETSC_ARCH=arch-debug<br class="">[0]PETSC ERROR: #1 TSHistoryGetLocFromTime() line 134 in /home/fbrarda/petsc/src/ts/interface/tshistory.c<br class="">[0]PETSC ERROR: #2 TSTrajectoryReconstruct_Private() line 55 in /home/fbrarda/petsc/src/ts/trajectory/utils/reconstruct.c<br class="">[0]PETSC ERROR: #3 TSTrajectoryGetVecs() line 239 in /home/fbrarda/petsc/src/ts/trajectory/interface/traj.c<br class=""><br class="">Is there any specific routine I can use to overcome this issue? Should I use VecScatter?<br class=""><br class="">I hope I made myself clear.<br class="">Best,<br class="">Francesco<br class=""></blockquote><br class=""></blockquote><br class=""></blockquote><br class=""></blockquote><br class=""></span></div></div></div></blockquote></div><br class=""></div></body></html>