<div dir="ltr"><div>Using DMDA to get automatic parallelization of a 0-D code is nonsense.</div><div>Ghost points are meant for spatially distributed data, not for 0D.<br></div><div>If you really want to use DM, you should use DMRedundant and call DMGlobalToLocal to go to your distributed to local full representation (it calls  MPI_Bcast internally).<br></div><div>My advice is to read the many examples we have.<br></div></div><br><div class="gmail_quote"><div dir="ltr" class="gmail_attr">Il giorno mar 4 mag 2021 alle ore 10:54 Francesco Brarda <<a href="mailto:brardafrancesco@gmail.com">brardafrancesco@gmail.com</a>> ha scritto:<br></div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div style="overflow-wrap: break-word;">Thank you very much everyone. I do have one more question. <br><blockquote type="cite">For what you want to do you can use 3,1,2.  This says three "spatial" points, 1 dof at each "spatial" point and 2 ghost points. In your case "spatial" does not mean spatial in space it is just three abstract points. <br></blockquote>In this case, since I have 2 ghost points, should I change the DMBoundaryType? <br>For one process this works, but with 2 or 3 procs it doesn’t. The error I have is the following:<div><br><font face="Menlo">Solving a non-linear TS problem, number of processors = 2<br>[1]PETSC ERROR: --------------------- Error Message --------------------------------------------------------------<br>[1]PETSC ERROR: Argument out of range<br>[1]PETSC ERROR: Local x-width of domain x 1 is smaller than stencil width s 2<br>[1]PETSC ERROR: See <a href="https://www.mcs.anl.gov/petsc/documentation/faq.html" target="_blank">https://www.mcs.anl.gov/petsc/documentation/faq.html</a> for trouble shooting.<br>[1]PETSC ERROR: Petsc Release Version 3.14.4, unknown <br>[1]PETSC ERROR: ./test_ic on a arch-debug named srvulx13 by fbrarda Mon May  3 17:46:37 2021<br>[1]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>[1]PETSC ERROR: #1 DMSetUp_DA_1D() line 199 in /home/fbrarda/petsc/src/dm/impls/da/da1.c<br>[1]PETSC ERROR: #2 DMSetUp_DA() line 20 in /home/fbrarda/petsc/src/dm/impls/da/dareg.c<br>[1]PETSC ERROR: #3 DMSetUp() line 787 in /home/fbrarda/petsc/src/dm/interface/dm.c<br>[1]PETSC ERROR: #4 main() line 232 in test_ic.c<br>[1]PETSC ERROR: No PETSc Option Table entries<br>[1]PETSC ERROR: ----------------End of Error Message -------send entire error message to <a href="mailto:petsc-maint@mcs.anl.gov" target="_blank">petsc-maint@mcs.anl.gov</a>----------<br></font><br><blockquote type="cite">Il giorno 2 mag 2021, alle ore 18:54, Barry Smith <<a href="mailto:bsmith@petsc.dev" target="_blank">bsmith@petsc.dev</a>> ha scritto:<br><br><br>[0]PETSC ERROR: Wrong subtype object:Parameter # 1 must have implementation da it is shell<br><br>Are you calling TSSetDM() to supply your created DMDA to the TS? Based on the error message you are not, it is using a default shell DM, which is what TS does if you do not provide them with a DM. You need to call TSSetDM() after you create the TS and the DMDA.<br><br><br><blockquote type="cite">On Apr 29, 2021, at 2:57 AM, Francesco Brarda <<a href="mailto:brardafrancesco@gmail.com" target="_blank">brardafrancesco@gmail.com</a>> wrote:<br><br>I defined the DM as follows<br>ierr = DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,3,3,3,NULL,&da);CHKERRQ(ierr);<br></blockquote><br>If you truly want one "spatial" point then you would want to use 1,3,0. This says 1 location in space, three degrees of freedom at that point and 0 ghost values (since there is only one spatial point there can be no ghost spatial values).  <br><br>BUT DMDA ALWAYS puts all degree of freedom at a point on the same process so this will not give you parallelism. All 3 dof will be on the same MPI rank.<br><br>For what you want to do you can use 3,1,2.  This says three "spatial" points, 1 dof at each "spatial" point and 2 ghost points. In your case "spatial" does not mean spatial in space it is just three abstract points. <br><br>The global vectors (from DMCreate or Get GlobalVector) will have 1 value on each rank. The local vector from DMCreate or Get LocalVector) will have three values on each rank. Your initial conditions code would be something like <br><br>if (rank == 0) {<br><blockquote type="cite">x[0] = appctx->N - appctx->p[2];<br></blockquote>} else if (rank == 1) {<br><blockquote type="cite">  x[1] = appctx->p[2];<br></blockquote>} else {<br><blockquote type="cite">  x[2] = 0.0;<br></blockquote>}<br><br><br>Your TSSetRHSFunction() would make a call to DMGetLocalVector(...&localX), do a DMGlobalToLocalBegin/End() from the input global X to localX, you would call DMDAVecGetArray(...,&xarray) on the localX and access all three values in xarray. The resulting computation of f the output vector would be something like <br><br>if (rank == 0) {<br><blockquote type="cite">farray[0] = your code that can use xarray[0], xarray[1], xarray[2]<br></blockquote>} else if (rank == 1) {<br><blockquote type="cite">farray[1] = your code that can use xarray[0], xarray[1], xarray[2]<br></blockquote>} else {<br><blockquote type="cite">farray[2] = your code that can use xarray[0], xarray[1], xarray[2]<br></blockquote>}<br><br>There are many examples of this pattern in the example tutorials.<br><br>When you implement a code with a spatial distribution you would use a dof of 3 at each point and not parallelize over the dof at each point. Likely you want to use DMNETWORK to manage the spatial distribution since it has a simple API and allows any number of different number of neighbors for each point. DMDA would not make sense  for true spatial distribution except in some truly trivial neighbor configurations. <br><br>Barry<br><br><br><br><blockquote type="cite"><br>I am not sure whether I correctly understood this command properly. The vector should have 3 components (S, I, R) and 3 DOF as it is defined only when the three coordinates have been set.<br>Then I create a global vector X. When I set the initial conditions as below <br><br>static PetscErrorCode InitialConditions(TS ts,Vec X, void *ctx)<br>{ <br>  PetscErrorCode    ierr;<br>  AppCtx            *appctx = (AppCtx*) ctx;<br>  PetscScalar       *x;<br>  DM                da;<br>  <br>  PetscFunctionBeginUser;<br>  ierr = TSGetDM(ts,&da);CHKERRQ(ierr);<br>  <br>  /* Get pointers to vector data */<br>  ierr = DMDAVecGetArray(da,X,(void*)&x);CHKERRQ(ierr);<br>  <br>  x[0] = appctx->N - appctx->p[2];<br>  x[1] = appctx->p[2];<br>  x[2] = 0.0;<br>  <br>  ierr = DMDAVecRestoreArray(da,X,(void*)&x);CHKERRQ(ierr);<br>  PetscFunctionReturn(0);<br>}<br><br>I have the error:<br><br>[0]PETSC ERROR: Invalid argument<br>[0]PETSC ERROR: Wrong subtype object:Parameter # 1 must have implementation da it is shell<br>[0]PETSC ERROR: See <a href="https://www.mcs.anl.gov/petsc/documentation/faq.html" target="_blank">https://www.mcs.anl.gov/petsc/documentation/faq.html</a> for trouble shooting.<br>[0]PETSC ERROR: Petsc Release Version 3.14.4, unknown <br>[0]PETSC ERROR: ./par_sir_model on a arch-debug named srvulx13 by fbrarda Thu Apr 29 09:36:17 2021<br>[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>[0]PETSC ERROR: #1 DMDAVecGetArray() line 48 in /home/fbrarda/petsc/src/dm/impls/da/dagetarray.c<br>[0]PETSC ERROR: #2 InitialConditions() line 175 in par_sir_model.c<br>[0]PETSC ERROR: #3 main() line 295 in par_sir_model.c<br>[0]PETSC ERROR: No PETSc Option Table entries<br><br>I would be very happy to receive any advices to fix the code.<br>Best,<br>Francesco<br><br><blockquote type="cite">Il giorno 20 apr 2021, alle ore 21:35, Matthew Knepley <<a href="mailto:knepley@gmail.com" target="_blank">knepley@gmail.com</a>> ha scritto:<br><br>On Tue, Apr 20, 2021 at 1:17 PM Francesco Brarda <<a href="mailto:brardafrancesco@gmail.com" target="_blank">brardafrancesco@gmail.com</a>> wrote:<br>Thank you for the advices, I would just like to convert the code I already have to see what might happen once parallelized.<br>Do you think it is better to put the 3 equations into a 1d Distributed Array with 3 dofs and run the job with multiple procs regardless of how many equations I have? Is it possible? <br><br>If you plan in the end to use a structured grid, this is a great plan. If not, this is not a good plan.<br><br>  Thanks,<br><br>     Matt<br> <br>Thank you,<br>Francesco<br><br><blockquote type="cite">Il giorno 20 apr 2021, alle ore 17:57, Stefano Zampini <<a href="mailto:stefano.zampini@gmail.com" target="_blank">stefano.zampini@gmail.com</a>> ha scritto:<br><br>It does not make sense to parallelize to 1 equation per process, unless that single equation per process is super super super costly.<br>Is this work you are doing used to understand PETSc parallelization strategy? if so, there are multiple examples in the sourcetree that you can look at to populate matrices and vectors in parallel<br><br>Il giorno mar 20 apr 2021 alle ore 17:52 Francesco Brarda <<a href="mailto:brardafrancesco@gmail.com" target="_blank">brardafrancesco@gmail.com</a>> ha scritto:<br>In principle the entire code was for 1 proc only. The functions were built with VecGetArray(). While adapting the code for multiple procs I thought using VecGetOwnershipRange was a possible way to allocate the equations in the vector using multiple procs. What do you think, please? <br><br>Thank you,<br>Francesco<br><br><blockquote type="cite">Il giorno 20 apr 2021, alle ore 16:43, Matthew Knepley <<a href="mailto:knepley@gmail.com" target="_blank">knepley@gmail.com</a>> ha scritto:<br><br>On Tue, Apr 20, 2021 at 10:41 AM Francesco Brarda <<a href="mailto:brardafrancesco@gmail.com" target="_blank">brardafrancesco@gmail.com</a>> wrote:<br>I was trying to follow Barry's advice some time ago, but I guess that's not the way he meant it. How should I refer to the values contained in x? With Distributed Arrays?<br><br>That is how you get values from x. However, I cannot understand at all what you are doing with "mybase".<br><br>   Matt<br> <br>Thanks <br>Francesco<br><br><blockquote type="cite"><blockquote type="cite"> Even though it will not scale and will deliver slower performance it is completely possible for you to solve the 3 variable problem using 3 MPI ranks. Or 10 mpi ranks. You would just create vectors/matrices with 1 degree of freedom for the first three ranks and no degrees of freedom for the later ranks. During your function evaluation (and Jacobian evaluation) for TS you will need to set up the appropriate communication to get the values you need on each rank to evaluate the parts of the function evaluation needed by that rank. This is true for parallelizing any computation.<br><br> Barry<br></blockquote></blockquote><br><br><br><blockquote type="cite">Il giorno 20 apr 2021, alle ore 15:40, Matthew Knepley <<a href="mailto:knepley@gmail.com" target="_blank">knepley@gmail.com</a>> ha scritto:<br><br>On Tue, Apr 20, 2021 at 9:36 AM Francesco Brarda <<a href="mailto:brardafrancesco@gmail.com" target="_blank">brardafrancesco@gmail.com</a>> wrote:<br>Hi!<br>I tried to implement the SIR model taking into account the fact that I will only use 3 MPI ranks at this moment.<br>I built vectors and matrices following the examples already available. In particular, I defined the functions required similarly (RHSFunction, IFunction, IJacobian), as follows:<br><br>I don't think this makes sense. You use "mybase" to distinguish between 3 procs, which would indicate that each procs has only<br>1 degree of freedom. However, you use x[1] on each proc, indicating it has at least 2 dofs.<br><br>  Thanks,<br><br>     Matt<br> <br>static PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ctx)<br>{ <br>  PetscErrorCode    ierr;<br>  AppCtx            *appctx = (AppCtx*) ctx;<br>  PetscScalar       f;//, *x_localptr; <br>  const PetscScalar *x;<br>  PetscInt          mybase;<br>  <br>  PetscFunctionBeginUser;<br>  ierr = VecGetOwnershipRange(X,&mybase,NULL);CHKERRQ(ierr);<br>  ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);<br>  if (mybase == 0) {<br>    f    = (PetscScalar) (-appctx->p1*x[0]*x[1]/appctx->N);<br>    ierr = VecSetValues(F,1,&mybase,&f,INSERT_VALUES);<br>  }<br>  if (mybase == 1) {<br>    f    = (PetscScalar) (appctx->p1*x[0]*x[1]/appctx->N-appctx->p2*x[1]);<br>    ierr = VecSetValues(F,1,&mybase,&f,INSERT_VALUES);<br>  }<br>  if (mybase == 2) {<br>    f    = (PetscScalar) (appctx->p2*x[1]);<br>    ierr = VecSetValues(F,1,&mybase,&f,INSERT_VALUES);<br>  }<br>  ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);<br>  ierr = VecAssemblyBegin(F);CHKERRQ(ierr);<br>  ierr = VecAssemblyEnd(F);CHKERRQ(ierr);<br>  PetscFunctionReturn(0);<br>}<br><br><br>Whilst for the Jacobian I did:<br><br><br>static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,void *ctx)<br>{ <br>  PetscErrorCode    ierr;<br>  AppCtx            *appctx = (AppCtx*) ctx;<br>  PetscInt          mybase, rowcol[] = {0,1,2};<br>  const PetscScalar *x;<br>  <br>  PetscFunctionBeginUser;<br>  ierr = MatGetOwnershipRange(B,&mybase,NULL);CHKERRQ(ierr);<br>  ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);<br>  if (mybase == 0) {<br>    const PetscScalar J[] = {a + appctx->p1*x[1]/appctx->N, appctx->p1*x[0]/appctx->N, 0};<br>    ierr = MatSetValues(B,1,&mybase,3,rowcol,J,INSERT_VALUES);CHKERRQ(ierr);<br>  }<br>  if (mybase == 1) {<br>    const PetscScalar J[] = {- appctx->p1*x[1]/appctx->N, a - appctx->p1*x[0]/appctx->N + appctx->p2, 0};<br>    ierr = MatSetValues(B,1,&mybase,3,rowcol,J,INSERT_VALUES);CHKERRQ(ierr);<br>  }<br>  if (mybase == 2) {<br>    const PetscScalar J[] = {0, - appctx->p2, a};<br>    ierr = MatSetValues(B,1,&mybase,3,rowcol,J,INSERT_VALUES);CHKERRQ(ierr);<br>  }<br>  ierr    = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);<br>  <br>  ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);<br>  ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);<br>  if (A != B) {<br>    ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);<br>    ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);<br>  }<br>  PetscFunctionReturn(0);<br>}<br><br>This code does not provide the correct result, that is, the solution is the initial condition, either using implicit or explicit methods. Is the way I defined these objects wrong? How can I fix it? <br>I also tried to print the Jacobian with the following commands but it does not work (blank rows and error message). How should I print the Jacobian?<br><br>ierr = TSGetIJacobian(ts,NULL,&K, NULL, NULL); CHKERRQ(ierr);<br>ierr = MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);<br>ierr = MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);<br>ierr = MatView(K,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);                            <br><br>I would very much appreciate any kind of help or advice.<br>Best,<br>Francesco<br><br><blockquote type="cite">Il giorno 2 apr 2021, alle ore 04:45, Barry Smith <<a href="mailto:bsmith@petsc.dev" target="_blank">bsmith@petsc.dev</a>> ha scritto:<br><br><br><br><blockquote type="cite" style="font-family:Helvetica;font-size:12px;font-style:normal;font-variant-caps:normal;font-weight:normal;letter-spacing:normal;text-align:start;text-indent:0px;text-transform:none;white-space:normal;word-spacing:0px;text-decoration:none">On Apr 1, 2021, at 9:17 PM, Zhang, Hong via petsc-users <<a href="mailto:petsc-users@mcs.anl.gov" target="_blank">petsc-users@mcs.anl.gov</a>> wrote:<br><br><br><br><blockquote type="cite">On Mar 31, 2021, at 2:53 AM, Francesco Brarda <<a href="mailto:brardafrancesco@gmail.com" target="_blank">brardafrancesco@gmail.com</a>> wrote:<br><br>Hi everyone!<br><br>I am trying to solve a system of 3 ODEs (a basic SIR model) with TS. Sequentially works pretty well, but I need to switch it into a parallel version. <br>I started working with TS not very long time ago, there are few questions I’d like to share with you and if you have any advices I’d be happy to hear.<br>First of all, do I need to use a DM object even if the model is only time dependent? All the examples I found were using that object for the other variable when solving PDEs.<br></blockquote><br>Are you considering SIR on a spatial domain? If so, you can parallelize your model in the spatial domain using DM. Splitting the three variables in the ODE among processors would not scale.<br></blockquote><br> Even though it will not scale and will deliver slower performance it is completely possible for you to solve the 3 variable problem using 3 MPI ranks. Or 10 mpi ranks. You would just create vectors/matrices with 1 degree of freedom for the first three ranks and no degrees of freedom for the later ranks. During your function evaluation (and Jacobian evaluation) for TS you will need to set up the appropriate communication to get the values you need on each rank to evaluate the parts of the function evaluation needed by that rank. This is true for parallelizing any computation.<br><br> Barry<br><br><br><br><br><blockquote type="cite" style="font-family:Helvetica;font-size:12px;font-style:normal;font-variant-caps:normal;font-weight:normal;letter-spacing:normal;text-align:start;text-indent:0px;text-transform:none;white-space:normal;word-spacing:0px;text-decoration:none"><br>Hong (Mr.)  <br><br><blockquote type="cite">When I preallocate the space for the Jacobian matrix, is it better to decide the local or global space?<br><br>Best,<br>Francesco<br></blockquote></blockquote></blockquote><br><br><br>-- <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><br><a href="https://www.cse.buffalo.edu/~knepley/" target="_blank">https://www.cse.buffalo.edu/~knepley/</a><br></blockquote><br><br><br>-- <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><br><a href="https://www.cse.buffalo.edu/~knepley/" target="_blank">https://www.cse.buffalo.edu/~knepley/</a><br></blockquote><br><br><br>-- <br>Stefano<br></blockquote><br><br><br>-- <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><br><a href="https://www.cse.buffalo.edu/~knepley/" target="_blank">https://www.cse.buffalo.edu/~knepley/</a><br></blockquote><br></blockquote><br></blockquote><br></div></div></blockquote></div><br clear="all"><br>-- <br><div dir="ltr" class="gmail_signature">Stefano</div>