<div dir="ltr"><br><br><div class="gmail_quote"><div dir="ltr">On Wed, Jul 4, 2018 at 1:09 PM Smith, Barry F. <<a href="mailto:bsmith@mcs.anl.gov">bsmith@mcs.anl.gov</a>> wrote:<br></div><blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex"><br>
<br>
> On Jul 4, 2018, at 9:40 AM, Mark Adams <<a href="mailto:mfadams@lbl.gov" target="_blank">mfadams@lbl.gov</a>> wrote:<br>
> <br>
> "-mat_seqaij_type seqaijmkl" just worked.<br>
<br>
    Mark,<br>
<br>
      Please clarify. Does this mean you can use -mat_seqaij_type seqaijmkl to satisfy all your needs now without changing any code?</blockquote><div><br></div><div>Apparently yes. I timed the code and it took much longer with NUM_OMP_THREADS > 1 (not sure what that is about but not my problem now), so that tells me that threads were activated.</div><div><br></div><div>I want to package up this simple approach and my branch that works with <span style="background-color:rgb(255,255,255);text-decoration-style:initial;text-decoration-color:initial;float:none;display:inline">-mat_type aijmkl<span> , and let the user analyze each. (So I am not worrying about performance, I just want correctness).</span></span></div><div> </div><blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex"> <br>
<br>
   Barry<br>
<br>
<br>
> On Wed, Jul 4, 2018 at 9:44 AM Mark Adams <<a href="mailto:mfadams@lbl.gov" target="_blank">mfadams@lbl.gov</a>> wrote:<br>
> <br>
> <br>
> On Wed, Jul 4, 2018 at 3:01 AM Richard Tran Mills <<a href="mailto:rtmills@anl.gov" target="_blank">rtmills@anl.gov</a>> wrote:<br>
> Hi Mark,<br>
> <br>
> I'd never looked at the code for MatCreateMPIAIJWithSeqAIJ(), but it looks like if you are using MPIAIJ matrices but you've set "-mat_seqaij_type seqaijmkl",<br>
> <br>
> No, I am using -mat_type aijmkl<br>
>  <br>
> then the "diagonal" portion of the MPIAIJ matrix ends up as a SEQAIJMKL instance, but the off-diagonal portion ends up as a plain SEQAIJ when it is created with MatCreateSeqAIJWithArrays(). (I have not stepped through this with a debugger to verify that this is actually what happens.) So, if we want the off-diagonal matrix to also be a SEQAIJMKL, then we either need to add some logic after the MatCreateSeqAIJWithArrays() call to check to see if a different subtype of SEQAIJ is being used and then convert it to the correct type, or we need a variant of MatCreateSeqAIJWithArrays() that will honor the subtype of SEQAIJ that we are using. Barry, since you added the concept of subtyping for SEQAIJ, what do you prefer?<br>
> <br>
> No preference.<br>
>  <br>
> I note that in some (most?) cases, I think we are OK or even better off with the off-diagonal matrix in the MPIAIJ not being a SEQAIJMKL one, because the PETSc routines can used the "compressed row" stuff to speed up handling the many zero rows that often end up in the off-diagonal, whereas the MKL routines cannot.<br>
> <br>
> <br>
> OK , I'll keep that in mind. I pushed a branch with fixes to get my gamg tests working (ksp ex56 primarily). I'm going to make a new branch and try "-mat_seqaij_type seqaijmkl".<br>
> <br>
> Thanks,<br>
> Mark<br>
>  <br>
> --Richard<br>
> <br>
> On Tue, Jul 3, 2018 at 4:31 PM, Mark Adams <<a href="mailto:mfadams@lbl.gov" target="_blank">mfadams@lbl.gov</a>> wrote:<br>
> I have having to fix AIJ methods that don't get the type from the args to set created matrix type, so as to keep the MKL typing.<br>
> <br>
> I am not sure how to fix this because I don't know to say 'make MPI from SEQ' in this method<br>
> <br>
> PetscErrorCode MatCreateMPIAIJWithSeqAIJ(MPI_Comm comm,Mat A,Mat B,const PetscInt garray[],Mat *mat)<br>
> {<br>
>  ....  <br>
>   ierr = MatCreate(comm,mat);CHKERRQ(ierr);<br>
>   ierr = MatGetSize(A,&m,&n);CHKERRQ(ierr);<br>
> .... <br>
>   ierr = MatSetSizes(*mat,m,n,PETSC_DECIDE,N);CHKERRQ(ierr);<br>
>   ierr = MatSetType(*mat,MATMPIAIJ);CHKERRQ(ierr);<br>
> <br>
> How should I do this?  Other places in the code I get the type from input (MPI) args and set created matrix types accordingly. I could just put a switch with an error for any new types that might come along in the future ...<br>
> <br>
> Mark<br>
> <br>
<br>
</blockquote></div></div>