<div dir="ltr">Ok thanks for clarifying that, i wasn't sure if there were different types,<div><br></div><div>Here is a stripped down version of my code, it seems like the preallocation is working now since the matrix population part is working without problem, but here it is for illustration purposes:</div><div><br></div><div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex">call DMSetMatrixPreallocateOnly(daDummy,PETSC_TRUE,ierr)<br>call DMCreateMatrix(daDummy,A,ierr)<br>call MatSetFromOptions(A,ierr)<br>call DMSetMatType(daDummy,MATMPIAIJVIENNACL,ierr)<br>call DMSetVecType(daDummy,VECMPIVIENNACL,ierr)<br>call MatMPIAIJSetPreallocation(A,19,PETSC_NULL_INTEGER,19,PETSC_NULL_INTEGER,ierr)<br>call MatSetUp(A,ierr)<br>[...]<br>            call MatSetValues(A,1,row,sumpos,pos(0:iter-1),vals(0:iter-1),INSERT_VALUES,ierr)<br>[...]<br>call MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY, ierr)<br>call MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY, ierr)</blockquote></div><div><br></div><div>Adding the first line there did the trick,</div><div><br></div><div>Now the problem seems to be the program is not recognizing the matrix as ViennaCL type when i try with more than one processor, i get now:</div><div><br></div><div>[0]PETSC ERROR: --------------------- Error Message --------------------------------------------------------------<br></div><div><div>[0]PETSC ERROR: No support for this operation for this object type</div><div>[0]PETSC ERROR: Currently only handles ViennaCL matrices</div><div>[0]PETSC ERROR: See <a href="http://www.mcs.anl.gov/petsc/documentation/faq.html">http://www.mcs.anl.gov/petsc/documentation/faq.html</a> for trouble shooting.</div><div>[0]PETSC ERROR: Petsc Development GIT revision: v3.9.2-549-g779ab53  GIT Date: 2018-05-31 17:31:13 +0300</div><div>[0]PETSC ERROR: ./gcmLEP.GPU on a cuda-debug named node50 by valera Wed Aug 15 14:44:22 2018</div><div>[0]PETSC ERROR: Configure options PETSC_ARCH=cuda-debug --with-mpi-dir=/usr/lib64/openmpi --COPTFLAGS=-O2 --CXXOPTFLAGS=-O2 --FOPTFLAGS=-O2 --with-shared-libraries=1 --with-debugging=1 --with-cuda=1 --CUDAFLAGS=-arch=sm_60 --with-blaslapack-dir=/usr/lib64 --download-viennacl</div><div>[0]PETSC ERROR: #1 PCSetUp_SAVIENNACL() line 47 in /home/valera/petsc/src/ksp/pc/impls/saviennaclcuda/<a href="http://saviennacl.cu">saviennacl.cu</a></div><div>[0]PETSC ERROR: #2 PCSetUp() line 932 in /home/valera/petsc/src/ksp/pc/interface/precon.c</div><div>[0]PETSC ERROR: #3 KSPSetUp() line 381 in /home/valera/petsc/src/ksp/ksp/interface/itfunc.c</div></div><div><br></div><div>When running with:<br></div><div><br></div><div>mpirun -n 1 ./gcmLEP.GPU tc=TestCases/LockRelease/LE_6x6x6/ jid=tiny_cuda_test_n2 -ksp_type cg -dm_vec_type viennacl -dm_mat_type aijviennacl -pc_type saviennacl -log_view</div><div><br></div><div><br></div><div>Thanks,</div><div><br></div><div><br></div><div><br></div><div><br></div><div><br></div><div><br></div><div><br></div><div><br></div><div><br></div></div><div class="gmail_extra"><br><div class="gmail_quote">On Wed, Aug 15, 2018 at 2:32 PM, Matthew Knepley <span dir="ltr"><<a href="mailto:knepley@gmail.com" target="_blank">knepley@gmail.com</a>></span> wrote:<br><blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex"><div dir="ltr"><div class="gmail_quote"><span class=""><div dir="ltr">On Wed, Aug 15, 2018 at 5:20 PM Manuel Valera <<a href="mailto:mvalera-w@sdsu.edu" target="_blank">mvalera-w@sdsu.edu</a>> wrote:<br></div><blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex"><div dir="ltr">It seems to be resumed on: I do not know how to preallocate a DM Matrix correctly.</div></blockquote><div><br></div></span><div>There is only one matrix type, Mat. There are no separate DM matrices. A DM can create a matrix for you</div><div>using DMCreateMatrix(), but that is a Mat and it is preallocated correctly. I am not sure what you are doing.</div><div><br></div><div>  Thanks,</div><div><br></div><div>    Matt</div><div><div class="h5"><div> </div><blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex"><div dir="ltr"><div>The interesting part is that it only breaks when i need to populate a GPU matrix from MPI, so kudos on that, but it seems i need to do better on my code to get this setup working,</div><div><br></div><div>Any help would be appreciated,</div><div><br></div><div>Thanks,</div><div><div><br></div><div><br></div></div></div><div class="gmail_extra"><br><div class="gmail_quote">On Wed, Aug 15, 2018 at 2:15 PM, Matthew Knepley <span dir="ltr"><<a href="mailto:knepley@gmail.com" target="_blank">knepley@gmail.com</a>></span> wrote:<br><blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex"><div dir="ltr"><div class="gmail_quote"><span><div dir="ltr">On Wed, Aug 15, 2018 at 4:53 PM Manuel Valera <<a href="mailto:mvalera-w@sdsu.edu" target="_blank">mvalera-w@sdsu.edu</a>> wrote:<br></div><blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex"><div dir="ltr">Thanks Matthew, <div><br></div><div>I try to do that when calling:</div><div><br></div><div><span style="font-size:12.8px;background-color:rgb(255,255,255);text-decoration-style:initial;text-decoration-color:initial;float:none;display:inline">call MatMPIAIJSetPreallocation(A,</span><span style="font-size:12.8px;background-color:rgb(255,255,255);text-decoration-style:initial;text-decoration-color:initial;float:none;display:inline">19<wbr>,PETSC_NULL_INTEGER,19,</span><span style="font-size:12.8px;background-color:rgb(255,255,255);text-decoration-style:initial;text-decoration-color:initial;float:none;display:inline">PETSC_<wbr>NULL_INTEGER,ierr)</span><br></div><div><span style="font-size:12.8px;background-color:rgb(255,255,255);text-decoration-style:initial;text-decoration-color:initial;float:none;display:inline"><br></span></div><div><span style="font-size:12.8px;background-color:rgb(255,255,255);text-decoration-style:initial;text-decoration-color:initial;float:none;display:inline">But i am not aware on how to do this for the DM if it needs something more specific/different,</span></div></div></blockquote><div><br></div></span><div>The error says that your preallocation is wrong for the values you are putting in. The DM does not control either,</div><div>so I do not understand your email.</div><div><br></div><div>  Thanks,</div><div><br></div><div>     Matt</div><div><div class="m_1335369059284202842m_6429250112161966131h5"><div> </div><blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex"><div dir="ltr"><div><span style="font-size:12.8px;background-color:rgb(255,255,255);text-decoration-style:initial;text-decoration-color:initial;float:none;display:inline">Thanks,</span></div></div><div class="gmail_extra"><br><div class="gmail_quote">On Wed, Aug 15, 2018 at 1:51 PM, Matthew Knepley <span dir="ltr"><<a href="mailto:knepley@gmail.com" target="_blank">knepley@gmail.com</a>></span> wrote:<br><blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex"><div dir="ltr"><div class="gmail_quote"><span><div dir="ltr">On Wed, Aug 15, 2018 at 4:39 PM Manuel Valera <<a href="mailto:mvalera-w@sdsu.edu" target="_blank">mvalera-w@sdsu.edu</a>> wrote:<br></div><blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex"><div dir="ltr">Hello PETSc devs,<div><br></div><div>I am running into an error when trying to use the MATMPIAIJVIENNACL Matrix type in MPI calls, the same code runs for MATSEQAIJVIENNACL type in one processor. The error happens when calling MatSetValues for this specific configuration. It does not occur when using MPI DMMatrix types only.</div></div></blockquote><div><br></div></span><div>The DM properly preallocates the matrix. I am assuming you do not here.</div><div><br></div><div>   Matt</div><div><div class="m_1335369059284202842m_6429250112161966131m_-5527426201967937475m_3347897161163976995h5"><div> </div><blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex"><div dir="ltr"><div>Any help will be appreciated, </div><div><br></div><div>Thanks,</div><div><br></div><div><br></div><div><br></div><div>My program call:</div><div><br></div><div>mpirun -n 2 ./gcmLEP.GPU tc=TestCases/LockRelease/LE_<wbr>6x6x6/ jid=tiny_cuda_test_n2 -ksp_type cg -dm_vec_type viennacl -dm_mat_type aijviennacl -pc_type saviennacl -log_view <br></div><div><br></div><div><br></div><div>The error (repeats after each MatSetValues call):</div><div><br></div><div><div>[1]PETSC ERROR: --------------------- Error Message ------------------------------<wbr>------------------------------<wbr>--</div><div>[1]PETSC ERROR: Argument out of range</div><div>[1]PETSC ERROR: Inserting a new nonzero at global row/column (75, 50) into matrix</div><div>[1]PETSC ERROR: See <a href="http://www.mcs.anl.gov/petsc/documentation/faq.html" target="_blank">http://www.mcs.anl.gov/petsc/<wbr>documentation/faq.html</a> for trouble shooting.</div><div>[1]PETSC ERROR: Petsc Development GIT revision: v3.9.2-549-g779ab53  GIT Date: 2018-05-31 17:31:13 +0300</div><div>[1]PETSC ERROR: ./gcmLEP.GPU on a cuda-debug named node50 by valera Wed Aug 15 13:10:44 2018</div><div>[1]PETSC ERROR: Configure options PETSC_ARCH=cuda-debug --with-mpi-dir=/usr/lib64/<wbr>openmpi --COPTFLAGS=-O2 --CXXOPTFLAGS=-O2 --FOPTFLAGS=-O2 --with-shared-libraries=1 --with-debugging=1 --with-cuda=1 --CUDAFLAGS=-arch=sm_60 --with-blaslapack-dir=/usr/<wbr>lib64 --download-viennacl</div><div>[1]PETSC ERROR: #1 MatSetValues_MPIAIJ() line 608 in /home/valera/petsc/src/mat/<wbr>impls/aij/mpi/mpiaij.c</div><div>[1]PETSC ERROR: #2 MatSetValues() line 1339 in /home/valera/petsc/src/mat/<wbr>interface/matrix.c</div></div><div><br></div><div><br></div><div>My Code structure:</div><div><br></div><div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex">call DMCreateMatrix(daDummy,A,ierr)<br>call MatSetFromOptions(A,ierr)<br>call MPI_Comm_size(PETSC_COMM_<wbr>WORLD, numprocs, ierr)<br>if (numprocs > 1) then  ! set matrix type parallel<br>    ! Get local size<br>    call DMDACreateNaturalVector(<wbr>daDummy,Tmpnat,ierr)<br>    call VecGetLocalSize(Tmpnat,<wbr>locsize,ierr)<br>    call VecDestroy(Tmpnat,ierr)<br>    ! Set matrix<br>#ifdef GPU<br>    call MatSetType(A,MATAIJVIENNACL,<wbr>ierr)<br>    call DMSetMatType(daDummy,<wbr>MATMPIAIJVIENNACL,ierr)<br>    call DMSetVecType(daDummy,<wbr>VECMPIVIENNACL,ierr)<br>    print*,'SETTING GPU TYPES'<br>#else<br>    call DMSetMatType(daDummy,<wbr>MATMPIAIJ,ierr)<br>    call DMSetMatType(daDummy,VECMPI,<wbr>ierr)<br>    call MatSetType(A,MATMPIAIJ,ierr)!<br>#endif<br>    call MatMPIAIJSetPreallocation(A,<wbr>19,PETSC_NULL_INTEGER,19,<wbr>PETSC_NULL_INTEGER,ierr)<br>else                    ! set matrix type sequential<br>#ifdef GPU<br>    call DMSetMatType(daDummy,<wbr>MATSEQAIJVIENNACL,ierr)<br>    call DMSetVecType(daDummy,<wbr>VECSEQVIENNACL,ierr)<br>    call MatSetType(A,<wbr>MATSEQAIJVIENNACL,ierr)<br>    print*,'SETTING GPU TYPES'<br>#else<br>    call DMSetMatType(daDummy,<wbr>MATSEQAIJ,ierr)<br>    call DMSetMatType(daDummy,VECSEQ,<wbr>ierr)<br>    call MatSetType(A,MATSEQAIJ,ierr)<br>#endif<br>call MatSetUp(A,ierr)<br>call getCenterInfo(daGrid,xstart,<wbr>ystart,zstart,xend,yend,zend)<br></blockquote><div> </div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex">do k=zstart,zend-1<br>    do j=ystart,yend-1<br>        do i=xstart,xend-1<br>[..]<br>           call MatSetValues(A,1,row,sumpos,<wbr>pos(0:iter-1),vals(0:iter-1),<wbr>INSERT_VALUES,ierr)<br>[..]</blockquote></div><div><br></div><div><br></div><div><br></div><div><br></div></div>
</blockquote></div></div></div><span class="m_1335369059284202842m_6429250112161966131m_-5527426201967937475m_3347897161163976995HOEnZb"><font color="#888888"><br clear="all"><div><br></div>-- <br><div dir="ltr" class="m_1335369059284202842m_6429250112161966131m_-5527426201967937475m_3347897161163976995m_3778373861107327029gmail_signature" data-smartmail="gmail_signature"><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.caam.rice.edu/~mk51/" target="_blank">https://www.cse.buffalo.edu/~<wbr>knepley/</a><br></div></div></div></div></div></font></span></div>
</blockquote></div><br></div>
</blockquote></div></div></div><div><div class="m_1335369059284202842m_6429250112161966131h5"><br clear="all"><div><br></div>-- <br><div dir="ltr" class="m_1335369059284202842m_6429250112161966131m_-5527426201967937475gmail_signature" data-smartmail="gmail_signature"><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.caam.rice.edu/~mk51/" target="_blank">https://www.cse.buffalo.edu/~<wbr>knepley/</a><br></div></div></div></div></div></div></div></div>
</blockquote></div><br></div>
</blockquote></div></div></div><div><div class="h5"><br clear="all"><div><br></div>-- <br><div dir="ltr" class="m_1335369059284202842gmail_signature" data-smartmail="gmail_signature"><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.caam.rice.edu/~mk51/" target="_blank">https://www.cse.buffalo.edu/~<wbr>knepley/</a><br></div></div></div></div></div></div></div></div>
</blockquote></div><br></div>