<html><head><meta http-equiv="content-type" content="text/html; charset=utf-8"></head><body style="overflow-wrap: break-word; -webkit-nbsp-mode: space; line-break: after-white-space;"><div><br></div>   I will add the two interfaces you requested today. (Likely you may need more also).<div><br></div><div>  Barry</div><div><br><div><br><blockquote type="cite"><div>On May 4, 2023, at 6:01 PM, Matthew Knepley <knepley@gmail.com> wrote:</div><br class="Apple-interchange-newline"><div><meta charset="UTF-8"><div dir="ltr" style="caret-color: rgb(0, 0, 0); font-family: Helvetica; font-size: 18px; font-style: normal; font-variant-caps: normal; font-weight: 400; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; -webkit-text-stroke-width: 0px; text-decoration: none;"><div dir="ltr">On Thu, May 4, 2023 at 1:43 PM LEONARDO MUTTI <<a href="mailto:leonardo.mutti01@universitadipavia.it" target="_blank">leonardo.mutti01@universitadipavia.it</a>> wrote:<br></div><div dir="ltr"><div class="gmail_quote"><blockquote class="gmail_quote" style="margin: 0px 0px 0px 0.8ex; border-left-width: 1px; border-left-style: solid; border-left-color: rgb(204, 204, 204); padding-left: 1ex;"><div dir="ltr"><div>Of course, I'll try to explain.</div><div><br></div><div>I am solving a parabolic equation with space-time FEM and I want an efficient solver/preconditioner for the resulting system.</div><div>The corresponding matrix, call it X, has an e.g. block bi-diagonal structure, if the cG(1)-dG(0) method is used (i.e. implicit Euler solved in batch).</div><div>Every block-row of X corresponds to a time instant.</div><div><br></div><div>I want to introduce parallelism in time by subdividing X into overlapping submatrices of e.g 2x2 or 3x3 blocks, along the block diagonal.</div><div>For instance, call X_i the individual blocks. The submatrices would be, for various i, (X_{i-1,i-1},X_{i-1,i};X_{i,i-1},X_{i,i}).</div><div>I'd like each submatrix to be solved in parallel, to combine the various results together in an ASM like fashion.</div><div>Every submatrix has thus a predecessor and a successor, and it overlaps with both, so that as far as I could understand, GASM has to be used in place of ASM.</div></div></blockquote><div><br></div><div>Yes, ordered that way you need GASM. I wonder if inverting the ordering would be useful, namely putting the time index on the inside.</div><div>Then the blocks would be over all time, but limited space, which is more the spirit of ASM I think.</div><div><br></div><div>Have you considered waveform relaxation for this problem?</div><div><br></div><div>   Thanks,</div><div><br></div><div>     Matt</div><div> <br></div><blockquote class="gmail_quote" style="margin: 0px 0px 0px 0.8ex; border-left-width: 1px; border-left-style: solid; border-left-color: rgb(204, 204, 204); padding-left: 1ex;"><div dir="ltr"><div>Hope this helps.</div><div>Best,<br></div><div>Leonardo</div><div><br></div><div class="gmail_quote"><div dir="ltr" class="gmail_attr">Il giorno gio 4 mag 2023 alle ore 18:05 Matthew Knepley <<a href="mailto:knepley@gmail.com" target="_blank">knepley@gmail.com</a>> ha scritto:<br></div><blockquote class="gmail_quote" style="margin: 0px 0px 0px 0.8ex; border-left-width: 1px; border-left-style: solid; border-left-color: rgb(204, 204, 204); padding-left: 1ex;"><div dir="ltr"><div dir="ltr">On Thu, May 4, 2023 at 11:24 AM LEONARDO MUTTI <<a href="mailto:leonardo.mutti01@universitadipavia.it" target="_blank">leonardo.mutti01@universitadipavia.it</a>> wrote:<br></div><div class="gmail_quote"><blockquote class="gmail_quote" style="margin: 0px 0px 0px 0.8ex; border-left-width: 1px; border-left-style: solid; border-left-color: rgb(204, 204, 204); padding-left: 1ex;"><div dir="ltr">Thank you for the help. <div>Adding to my example:<div><i>     <span class="Apple-converted-space"> </span>call PCGASMSetSubdomains(pc,NSub, subdomains_IS, inflated_IS,ierr)<br>     <span class="Apple-converted-space"> </span>call PCGASMDestroySubdomains(NSub,subdomains_IS,inflated_IS,ierr)<br></i></div><div>results in:<br></div><div><div><i>     <span class="Apple-converted-space"> </span>Error LNK2019 unresolved external symbol PCGASMDESTROYSUBDOMAINS referenced in function ...  <br></i></div><div><i>     <span class="Apple-converted-space"> </span>Error LNK2019 unresolved external symbol PCGASMSETSUBDOMAINS referenced in function ...  <br></i></div><div>I'm not sure if the interfaces are missing or if I have a compilation problem.</div></div></div></div></blockquote><div><br></div><div>I just want to make sure you really want GASM. It sounded like you might able to do what you want just with ASM.</div><div>Can you tell me again what you want to do overall?</div><div><br></div><div> <span class="Apple-converted-space"> </span>Thanks,</div><div><br></div><div>   <span class="Apple-converted-space"> </span>Matt</div><div> </div><blockquote class="gmail_quote" style="margin: 0px 0px 0px 0.8ex; border-left-width: 1px; border-left-style: solid; border-left-color: rgb(204, 204, 204); padding-left: 1ex;"><div dir="ltr"><div><div>Thank you again. </div><div>Best,<br>Leonardo</div></div></div><br><div class="gmail_quote"><div dir="ltr" class="gmail_attr">Il giorno sab 29 apr 2023 alle ore 20:30 Barry Smith <<a href="mailto:bsmith@petsc.dev" target="_blank">bsmith@petsc.dev</a>> ha scritto:<br></div><blockquote class="gmail_quote" style="margin: 0px 0px 0px 0.8ex; border-left-width: 1px; border-left-style: solid; border-left-color: rgb(204, 204, 204); padding-left: 1ex;"><div><div><br></div>   Thank you for the test code. I have a fix in the branch <span style="color: rgb(51, 50, 56); font-family: "GitLab Sans", -apple-system, BlinkMacSystemFont, "Segoe UI", Roboto, "Noto Sans", Ubuntu, Cantarell, "Helvetica Neue", sans-serif, "Apple Color Emoji", "Segoe UI Emoji", "Segoe UI Symbol", "Noto Color Emoji"; font-size: 14px; background-color: rgb(251, 250, 253);"> </span><a href="https://gitlab.com/petsc/petsc/-/commits/barry/2023-04-29/fix-pcasmcreatesubdomains2d" target="_blank" style="box-sizing: border-box; font-variant-ligatures: none; color: rgb(31, 117, 203); text-decoration: none; font-size: 13.3px;">barry/2023-04-29/fix-pcasmcreatesubdomains2d</a> with merge request <a href="https://gitlab.com/petsc/petsc/-/merge_requests/6394" target="_blank">https://gitlab.com/petsc/petsc/-/merge_requests/6394</a><div><br></div><div>   The functions did not have proper Fortran stubs and interfaces so I had to provide them manually in the new branch.<br><div><br></div><div>   Use</div><div><br></div><div>   git fetch</div><div>   git checkout <a href="https://gitlab.com/petsc/petsc/-/commits/barry/2023-04-29/fix-pcasmcreatesubdomains2d" target="_blank" style="box-sizing: border-box; font-variant-ligatures: none; color: rgb(31, 117, 203); text-decoration: none; font-size: 13.3px;">barry/2023-04-29/fix-pcasmcreatesubdomains2d</a></div><div>   ./configure etc</div><div><br></div><div>   Your now working test code is in src/ksp/ksp/tests/ex71f.F90  I had to change things slightly and I updated the error handling for the latest version.</div><div><br></div><div>   Please let us know if you have any later questions.</div><div><br></div><div> <span class="Apple-converted-space"> </span>Barry</div><div><br></div><div><br></div><div><br></div><div><div><br><blockquote type="cite"><div>On Apr 28, 2023, at 12:07 PM, LEONARDO MUTTI <<a href="mailto:leonardo.mutti01@universitadipavia.it" target="_blank">leonardo.mutti01@universitadipavia.it</a>> wrote:</div><br><div><div dir="ltr">Hello. I am having a hard time understanding the index sets to feed PCGASMSetSubdomains, and I am working in Fortran (as a PETSc novice). To get more intuition on how the IS objects behave I tried the following minimal (non) working example, which should tile a 16x16 matrix into 16 square, non-overlapping submatrices:<div><br></div><div>#include <petsc/finclude/petscmat.h><br></div><div>#include <petsc/finclude/petscksp.h><br>#include <petsc/finclude/petscpc.h><br>     <span class="Apple-converted-space"> </span>USE petscmat<br>     <span class="Apple-converted-space"> </span>USE petscksp<br>     <span class="Apple-converted-space"> </span>USE petscpc<br>      <br>     <span class="Apple-converted-space"> </span>Mat   :: A<br>     <span class="Apple-converted-space"> </span>PetscInt :: M, NSubx, dof, overlap, NSub<br>     <span class="Apple-converted-space"> </span>INTEGER :: I,J<br>     <span class="Apple-converted-space"> </span>PetscErrorCode      :: ierr<br>     <span class="Apple-converted-space"> </span>PetscScalar :: v<br>     <span class="Apple-converted-space"> </span>KSP            :: ksp<br>     <span class="Apple-converted-space"> </span>PC             :: pc<br>     <span class="Apple-converted-space"> </span>IS :: subdomains_IS, inflated_IS<br>     <span class="Apple-converted-space"> </span><br>     <span class="Apple-converted-space"> </span>call PetscInitialize(PETSC_NULL_CHARACTER , ierr)</div><div>      <br>!-----Create a dummy matrix<br>     <span class="Apple-converted-space"> </span>M = 16<br>     <span class="Apple-converted-space"> </span>call MatCreateAIJ(MPI_COMM_WORLD, PETSC_DECIDE, PETSC_DECIDE,<br>     & M, M,<br>     & PETSC_DEFAULT_INTEGER, PETSC_NULL_INTEGER,<br>     & PETSC_DEFAULT_INTEGER, PETSC_NULL_INTEGER,<br>     & A, ierr)<br>     <span class="Apple-converted-space"> </span><br>     <span class="Apple-converted-space"> </span>DO I=1,M<br>         DO J=1,M<br>           <span class="Apple-converted-space"> </span>v = I*J<br>           <span class="Apple-converted-space"> </span>CALL MatSetValue (A,I-1,J-1,v,<br>     &                     INSERT_VALUES , ierr)<br>         END DO<br>     <span class="Apple-converted-space"> </span>END DO<br>     <span class="Apple-converted-space"> </span><br>     <span class="Apple-converted-space"> </span>call MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY , ierr)<br>     <span class="Apple-converted-space"> </span>call MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY , ierr)<br>      </div><div>!-----Create KSP and PC<br>     <span class="Apple-converted-space"> </span>call KSPCreate(PETSC_COMM_WORLD,ksp, ierr)<br>     <span class="Apple-converted-space"> </span>call KSPSetOperators(ksp,A,A, ierr)<br>     <span class="Apple-converted-space"> </span>call KSPSetType(ksp,"bcgs",ierr)<br>     <span class="Apple-converted-space"> </span>call KSPGetPC(ksp,pc,ierr)<br>     <span class="Apple-converted-space"> </span>call KSPSetUp(ksp, ierr)<br>     <span class="Apple-converted-space"> </span>call PCSetType(pc,PCGASM, ierr)<br>     <span class="Apple-converted-space"> </span>call PCSetUp(pc , ierr)<br>     <span class="Apple-converted-space"> </span><br>!-----GASM setup<br>      NSubx = 4<br>     <span class="Apple-converted-space"> </span>dof = 1<br>     <span class="Apple-converted-space"> </span>overlap = 0<br>     <span class="Apple-converted-space"> </span><br>     <span class="Apple-converted-space"> </span>call PCGASMCreateSubdomains2D(pc,<span class="Apple-converted-space"> </span><br>     &      M, M,<br>     &      NSubx, NSubx,<br>     &      dof, overlap,<br>     &      NSub, subdomains_IS, inflated_IS, ierr)<br><br>     <span class="Apple-converted-space"> </span>call ISView(subdomains_IS, PETSC_VIEWER_STDOUT_WORLD, ierr)<br>     <span class="Apple-converted-space"> </span><br>     <span class="Apple-converted-space"> </span>call KSPDestroy(ksp, ierr)  <br>     <span class="Apple-converted-space"> </span>call PetscFinalize(ierr)</div><div><br></div><div>Running this on one processor, I get NSub = 4.</div><div>If PCASM and PCASMCreateSubdomains2D are used instead, I get NSub = 16 as expected.</div><div>Moreover, I get in the end "forrtl: severe (157): Program Exception - access violation". So:</div><div>1) why do I get two different results with ASM, and GASM?</div><div>2) why do I get access violation and how can I solve this?</div><div>In fact, in C, subdomains_IS, inflated_IS should pointers to IS objects. As I see on the Fortran interface, the arguments to PCGASMCreateSubdomains2D are IS objects:<br></div><div><br></div><div>       subroutine PCGASMCreateSubdomains2D(a,b,c,d,e,f,g,h,i,j,z)<br>       import tPC,tIS<br>       PC a ! PC<br>       PetscInt b ! PetscInt<br>       PetscInt c ! PetscInt<br>       PetscInt d ! PetscInt<br>       PetscInt e ! PetscInt<br>       PetscInt f ! PetscInt<br>       PetscInt g ! PetscInt<br>       PetscInt h ! PetscInt<br>       IS i ! IS<br>       IS j ! IS<br>       PetscErrorCode z<br>       end subroutine PCGASMCreateSubdomains2D<br></div><div>Thus:</div><div>3) what should be inside e.g., subdomains_IS? I expect it to contain, for every created subdomain, the list of rows and columns defining the subblock in the matrix, am I right?</div><div><br></div><div>Context: I have a block-tridiagonal system arising from space-time finite elements, and I want to solve it with GMRES+PCGASM preconditioner, where each overlapping submatrix is on the diagonal and of size 3x3 blocks (and spanning multiple processes). This is PETSc 3.17.1 on Windows.</div><div><br></div><div>Thanks in advance,</div><div>Leonardo</div></div></div></blockquote></div><br></div></div></div></blockquote></div></blockquote></div><br clear="all"><div><br></div><span>--<span class="Apple-converted-space"> </span></span><br><div dir="ltr"><div dir="ltr"><div><div dir="ltr"><div><div dir="ltr"><div>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</div><div><br></div><div><a href="http://www.cse.buffalo.edu/~knepley/" target="_blank">https://www.cse.buffalo.edu/~knepley/</a><br></div></div></div></div></div></div></div></div></blockquote></div></div></blockquote></div><br clear="all"><div><br></div><span>--<span class="Apple-converted-space"> </span></span><br><div dir="ltr"><div dir="ltr"><div><div dir="ltr"><div><div dir="ltr"><div>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</div><div><br></div><div><a href="http://www.cse.buffalo.edu/~knepley/" target="_blank">https://www.cse.buffalo.edu/~knepley/</a></div></div></div></div></div></div></div></div></div></div></blockquote></div><br></div></body></html>