<div dir="ltr"><div class="gmail_quote"><div dir="ltr">On Wed, Aug 15, 2018 at 4:39 PM Manuel Valera <<a href="mailto:mvalera-w@sdsu.edu">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><div>The DM properly preallocates the matrix. I am assuming you do not here.</div><div><br></div><div>   Matt</div><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_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 --------------------------------------------------------------</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/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/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>[1]PETSC ERROR: #1 MatSetValues_MPIAIJ() line 608 in /home/valera/petsc/src/mat/impls/aij/mpi/mpiaij.c</div><div>[1]PETSC ERROR: #2 MatSetValues() line 1339 in /home/valera/petsc/src/mat/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_WORLD, numprocs, ierr)<br>if (numprocs > 1) then  ! set matrix type parallel<br>    ! Get local size<br>    call DMDACreateNaturalVector(daDummy,Tmpnat,ierr)<br>    call VecGetLocalSize(Tmpnat,locsize,ierr)<br>    call VecDestroy(Tmpnat,ierr)<br>    ! Set matrix<br>#ifdef GPU<br>    call MatSetType(A,MATAIJVIENNACL,ierr)<br>    call DMSetMatType(daDummy,MATMPIAIJVIENNACL,ierr)<br>    call DMSetVecType(daDummy,VECMPIVIENNACL,ierr)<br>    print*,'SETTING GPU TYPES'<br>#else<br>    call DMSetMatType(daDummy,MATMPIAIJ,ierr)<br>    call DMSetMatType(daDummy,VECMPI,ierr)<br>    call MatSetType(A,MATMPIAIJ,ierr)!<br>#endif<br>    call MatMPIAIJSetPreallocation(A,19,PETSC_NULL_INTEGER,19,PETSC_NULL_INTEGER,ierr)<br>else                    ! set matrix type sequential<br>#ifdef GPU<br>    call DMSetMatType(daDummy,MATSEQAIJVIENNACL,ierr)<br>    call DMSetVecType(daDummy,VECSEQVIENNACL,ierr)<br>    call MatSetType(A,MATSEQAIJVIENNACL,ierr)<br>    print*,'SETTING GPU TYPES'<br>#else<br>    call DMSetMatType(daDummy,MATSEQAIJ,ierr)<br>    call DMSetMatType(daDummy,VECSEQ,ierr)<br>    call MatSetType(A,MATSEQAIJ,ierr)<br>#endif<br>call MatSetUp(A,ierr)<br>call getCenterInfo(daGrid,xstart,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,pos(0:iter-1),vals(0:iter-1),INSERT_VALUES,ierr)<br>[..]</blockquote></div><div><br></div><div><br></div><div><br></div><div><br></div></div>
</blockquote></div><br clear="all"><div><br></div>-- <br><div dir="ltr" class="gmail_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/~knepley/</a><br></div></div></div></div></div></div>