<div dir="ltr">Hello everyone,<div><br></div><div>I just had time to work on this again, and checked the code for errors on the matrix entries, this is the exact code i was using for creating the matrix without DMSetMatrixPreallocateOnly, using MatMPIAIJSetPreallocation and it worked that way, but trying this way i get the same 'Column too large' error using any number at the column position of MatSetValues, </div><div><br></div><div>I have set up my code to print the column argument (n) of MatSetValues and in this case is 7 (lower than 124), it still gives error, even entering a specific number in the MatSetValues column argument position gives the same error.</div><div><br></div><div>So next i went back to ex.50 here: <a href="http://www.mcs.anl.gov/petsc/petsc-current/src/ts/examples/tutorials/ex50.c.html">http://www.mcs.anl.gov/petsc/petsc-current/src/ts/examples/tutorials/ex50.c.html</a> and it has a very similar structure except the PetscMalloc1() call, so i tried adding that and got:</div><div><br></div><div> /home/valera/ParGCCOM/Src/DMDALaplacian.f90:114: undefined reference to `petscmalloc1_'</div><div><br></div><div>Any ideas on this behaviour?</div><div><br></div><div>Thanks so much,</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 Thu, Aug 16, 2018 at 11:20 AM, Smith, Barry F. <span dir="ltr"><<a href="mailto:bsmith@mcs.anl.gov" target="_blank">bsmith@mcs.anl.gov</a>></span> wrote:<br><blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex"><span class=""><br>
Column too large: col 10980 max 124<br>
<br>
</span>   You need to check the code that is generating the matrix entries. The matrix has 124 columns but you are attempting to put a value at column 10980<br>
<span class="HOEnZb"><font color="#888888"><br>
   Barry<br>
</font></span><div class="HOEnZb"><div class="h5"><br>
<br>
> On Aug 15, 2018, at 9:44 PM, Manuel Valera <<a href="mailto:mvalera-w@sdsu.edu">mvalera-w@sdsu.edu</a>> wrote:<br>
> <br>
> Thanks Matthew and Barry,<br>
> <br>
> Now my code looks like:<br>
> <br>
> call DMSetMatrixPreallocateOnly(<wbr>daDummy,PETSC_TRUE,ierr)<br>
> call DMSetMatType(daDummy,<wbr>MATMPIAIJVIENNACL,ierr)<br>
> call DMSetVecType(daDummy,<wbr>VECMPIVIENNACL,ierr)<br>
> call DMCreateMatrix(daDummy,A,ierr)<br>
> call MatSetFromOptions(A,ierr)<br>
> call MatSetUp(A,ierr)<br>
> [...]<br>
>             call MatSetValues(A,1,row,sumpos,<wbr>pos(0:iter-1),vals(0:iter-1),<wbr>INSERT_VALUES,ierr)<br>
> [...]<br>
> call MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY, ierr)<br>
> call MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY, ierr)<br>
> <br>
> And i get a different error, now is:<br>
> <br>
> [0]PETSC ERROR: --------------------- Error Message ------------------------------<wbr>------------------------------<wbr>--<br>
> [0]PETSC ERROR: Argument out of range<br>
> [0]PETSC ERROR: Column too large: col 10980 max 124<br>
> [0]PETSC ERROR: See <a href="http://www.mcs.anl.gov/petsc/documentation/faq.html" rel="noreferrer" target="_blank">http://www.mcs.anl.gov/petsc/<wbr>documentation/faq.html</a> for trouble shooting.<br>
> [0]PETSC ERROR: Petsc Development GIT revision: v3.9.2-549-g779ab53  GIT Date: 2018-05-31 17:31:13 +0300<br>
> [0]PETSC ERROR: ./gcmLEP.GPU on a cuda-debug named node50 by valera Wed Aug 15 19:40:00 2018<br>
> [0]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<br>
> [0]PETSC ERROR: #1 MatSetValues_SeqAIJ() line 442 in /home/valera/petsc/src/mat/<wbr>impls/aij/seq/aij.c<br>
> [0]PETSC ERROR: #2 MatSetValues() line 1339 in /home/valera/petsc/src/mat/<wbr>interface/matrix.c<br>
> <br>
> <br>
> Thanks again,<br>
> <br>
> <br>
> <br>
> <br>
> <br>
> <br>
> <br>
> <br>
> On Wed, Aug 15, 2018 at 7:02 PM, Smith, Barry F. <<a href="mailto:bsmith@mcs.anl.gov">bsmith@mcs.anl.gov</a>> wrote:<br>
> <br>
>   Should be<br>
> <br>
> call DMSetMatType(daDummy,<wbr>MATMPIAIJVIENNACL,ierr)<br>
> call DMSetVecType(daDummy,<wbr>VECMPIVIENNACL,ierr)<br>
> call DMCreateMatrix(daDummy,A,ierr)<br>
> <br>
>   and remove the rest. You need to set the type of Mat you want the DM to return BEFORE you create the matrix.<br>
> <br>
>   Barry<br>
> <br>
> <br>
> <br>
> > On Aug 15, 2018, at 4:45 PM, Manuel Valera <<a href="mailto:mvalera-w@sdsu.edu">mvalera-w@sdsu.edu</a>> wrote:<br>
> > <br>
> > Ok thanks for clarifying that, i wasn't sure if there were different types,<br>
> > <br>
> > 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:<br>
> > <br>
> > call DMSetMatrixPreallocateOnly(<wbr>daDummy,PETSC_TRUE,ierr)<br>
> > call DMCreateMatrix(daDummy,A,ierr)<br>
> > call MatSetFromOptions(A,ierr)<br>
> > call DMSetMatType(daDummy,<wbr>MATMPIAIJVIENNACL,ierr)<br>
> > call DMSetVecType(daDummy,<wbr>VECMPIVIENNACL,ierr)<br>
> > call MatMPIAIJSetPreallocation(A,<wbr>19,PETSC_NULL_INTEGER,19,<wbr>PETSC_NULL_INTEGER,ierr)<br>
> > call MatSetUp(A,ierr)<br>
> > [...]<br>
> >             call MatSetValues(A,1,row,sumpos,<wbr>pos(0:iter-1),vals(0:iter-1),<wbr>INSERT_VALUES,ierr)<br>
> > [...]<br>
> > call MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY, ierr)<br>
> > call MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY, ierr)<br>
> > <br>
> > Adding the first line there did the trick,<br>
> > <br>
> > 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:<br>
> > <br>
> > [0]PETSC ERROR: --------------------- Error Message ------------------------------<wbr>------------------------------<wbr>--<br>
> > [0]PETSC ERROR: No support for this operation for this object type<br>
> > [0]PETSC ERROR: Currently only handles ViennaCL matrices<br>
> > [0]PETSC ERROR: See <a href="http://www.mcs.anl.gov/petsc/documentation/faq.html" rel="noreferrer" target="_blank">http://www.mcs.anl.gov/petsc/<wbr>documentation/faq.html</a> for trouble shooting.<br>
> > [0]PETSC ERROR: Petsc Development GIT revision: v3.9.2-549-g779ab53  GIT Date: 2018-05-31 17:31:13 +0300<br>
> > [0]PETSC ERROR: ./gcmLEP.GPU on a cuda-debug named node50 by valera Wed Aug 15 14:44:22 2018<br>
> > [0]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<br>
> > [0]PETSC ERROR: #1 PCSetUp_SAVIENNACL() line 47 in /home/valera/petsc/src/ksp/pc/<wbr>impls/saviennaclcuda/<a href="http://saviennacl.cu" rel="noreferrer" target="_blank">saviennac<wbr>l.cu</a><br>
> > [0]PETSC ERROR: #2 PCSetUp() line 932 in /home/valera/petsc/src/ksp/pc/<wbr>interface/precon.c<br>
> > [0]PETSC ERROR: #3 KSPSetUp() line 381 in /home/valera/petsc/src/ksp/<wbr>ksp/interface/itfunc.c<br>
> > <br>
> > When running with:<br>
> > <br>
> > mpirun -n 1 ./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>
> > <br>
> > <br>
> > Thanks,<br>
> > <br>
> > <br>
> > <br>
> > <br>
> > <br>
> > <br>
> > <br>
> > <br>
> > <br>
> > <br>
> > On Wed, Aug 15, 2018 at 2:32 PM, Matthew Knepley <<a href="mailto:knepley@gmail.com">knepley@gmail.com</a>> wrote:<br>
> > On Wed, Aug 15, 2018 at 5:20 PM Manuel Valera <<a href="mailto:mvalera-w@sdsu.edu">mvalera-w@sdsu.edu</a>> wrote:<br>
> > It seems to be resumed on: I do not know how to preallocate a DM Matrix correctly.<br>
> > <br>
> > There is only one matrix type, Mat. There are no separate DM matrices. A DM can create a matrix for you<br>
> > using DMCreateMatrix(), but that is a Mat and it is preallocated correctly. I am not sure what you are doing.<br>
> > <br>
> >   Thanks,<br>
> > <br>
> >     Matt<br>
> >  <br>
> > 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,<br>
> > <br>
> > Any help would be appreciated,<br>
> > <br>
> > Thanks,<br>
> > <br>
> > <br>
> > <br>
> > On Wed, Aug 15, 2018 at 2:15 PM, Matthew Knepley <<a href="mailto:knepley@gmail.com">knepley@gmail.com</a>> wrote:<br>
> > On Wed, Aug 15, 2018 at 4:53 PM Manuel Valera <<a href="mailto:mvalera-w@sdsu.edu">mvalera-w@sdsu.edu</a>> wrote:<br>
> > Thanks Matthew, <br>
> > <br>
> > I try to do that when calling:<br>
> > <br>
> > call MatMPIAIJSetPreallocation(A,<wbr>19,PETSC_NULL_INTEGER,19,<wbr>PETSC_NULL_INTEGER,ierr)<br>
> > <br>
> > But i am not aware on how to do this for the DM if it needs something more specific/different,<br>
> > <br>
> > The error says that your preallocation is wrong for the values you are putting in. The DM does not control either,<br>
> > so I do not understand your email.<br>
> > <br>
> >   Thanks,<br>
> > <br>
> >      Matt<br>
> >  <br>
> > Thanks,<br>
> > <br>
> > On Wed, Aug 15, 2018 at 1:51 PM, Matthew Knepley <<a href="mailto:knepley@gmail.com">knepley@gmail.com</a>> wrote:<br>
> > 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>
> > Hello PETSc devs,<br>
> > <br>
> > 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.<br>
> > <br>
> > The DM properly preallocates the matrix. I am assuming you do not here.<br>
> > <br>
> >    Matt<br>
> >  <br>
> > Any help will be appreciated, <br>
> > <br>
> > Thanks,<br>
> > <br>
> > <br>
> > <br>
> > My program call:<br>
> > <br>
> > 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>
> > <br>
> > <br>
> > The error (repeats after each MatSetValues call):<br>
> > <br>
> > [1]PETSC ERROR: --------------------- Error Message ------------------------------<wbr>------------------------------<wbr>--<br>
> > [1]PETSC ERROR: Argument out of range<br>
> > [1]PETSC ERROR: Inserting a new nonzero at global row/column (75, 50) into matrix<br>
> > [1]PETSC ERROR: See <a href="http://www.mcs.anl.gov/petsc/documentation/faq.html" rel="noreferrer" target="_blank">http://www.mcs.anl.gov/petsc/<wbr>documentation/faq.html</a> for trouble shooting.<br>
> > [1]PETSC ERROR: Petsc Development GIT revision: v3.9.2-549-g779ab53  GIT Date: 2018-05-31 17:31:13 +0300<br>
> > [1]PETSC ERROR: ./gcmLEP.GPU on a cuda-debug named node50 by valera Wed Aug 15 13:10:44 2018<br>
> > [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<br>
> > [1]PETSC ERROR: #1 MatSetValues_MPIAIJ() line 608 in /home/valera/petsc/src/mat/<wbr>impls/aij/mpi/mpiaij.c<br>
> > [1]PETSC ERROR: #2 MatSetValues() line 1339 in /home/valera/petsc/src/mat/<wbr>interface/matrix.c<br>
> > <br>
> > <br>
> > My Code structure:<br>
> > <br>
> > 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>
> >  <br>
> > 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>
> > [..]<br>
> > <br>
> > <br>
> > <br>
> > <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/" rel="noreferrer" target="_blank">https://www.cse.buffalo.edu/~<wbr>knepley/</a><br>
> > <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/" rel="noreferrer" target="_blank">https://www.cse.buffalo.edu/~<wbr>knepley/</a><br>
> > <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/" rel="noreferrer" target="_blank">https://www.cse.buffalo.edu/~<wbr>knepley/</a><br>
> > <br>
> <br>
> <br>
<br>
</div></div></blockquote></div><br></div>