<div dir="ltr"><div class="gmail_extra"><div class="gmail_quote">On Tue, Nov 25, 2014 at 12:54 PM, Kirill Voronin <span dir="ltr"><<a href="mailto:kvoronin@labchem.sscc.ru" target="_blank">kvoronin@labchem.sscc.ru</a>></span> wrote:<br><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left-width:1px;border-left-color:rgb(204,204,204);border-left-style:solid;padding-left:1ex"><br>
Thank you for the prompt answer!<br>
<br>
I shifted the indices but when calling MatSeqSBAIJSetPreallocationCSR got<br>
the error - nnz cannot be greater than block row length: local row 0 value<br>
8 rowlength 1.<br></blockquote><div><br></div><div>It sounds like the input you are giving is for actual rows, but this call expects the</div><div>CSR to reference blocks. Does that make sense?</div><div> </div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left-width:1px;border-left-color:rgb(204,204,204);border-left-style:solid;padding-left:1ex">
Can anyone please make a bit more clear some details of block AIJ format<br>
in PETSC?<br>
<br>
My code looks like:<br>
-----------------------------------------------------<br>
ierr = MatCreate(comm,&Afromcsr);CHKERRQ(ierr);<br>
ierr = MatSetType(Afromcsr,"sbaij"); CHKERRQ(ierr);<br>
ierr = MatSetSizes(Afromcsr,dim,dim,dim,dim);CHKERRQ(ierr);<br>
ierr = MatSetFromOptions(Afromcsr);CHKERRQ(ierr);<br>
<br>
printf("fill in for values of acsr \n");<br>
<br>
for (i = 0; i <= dim; i++)<br>
ia[i] -= 1;<br>
for (i = 0; i < ia[dim]; i++)<br>
ja[i] -= 1;<br>
ierr = MatSeqSBAIJSetPreallocationCSR(Afromcsr, dim, ia, ja, acsr);<br>
CHKERRQ(ierr);<br>
-----------------------------------------------------------<br>
The questions:<br>
<br>
1) Why no MatMPISBAIJSetPreaalocationCSR? Or should<br>
MatCreateMPISBAIJWithArrays be used instead?<br></blockquote><div><br></div><div>Not written yet.</div><div> </div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left-width:1px;border-left-color:rgb(204,204,204);border-left-style:solid;padding-left:1ex">
2) What is block size and why it is assumed to be 1 (according to the<br>
error I got).<br></blockquote><div><br></div><div>Block size is the size of the dense blocks out of which the matrix is constructed.</div><div>This can achieve better performance since indexing is only done on blocks, not</div><div>individual elements, and dense algebra can be done on block-block computation.</div><div> </div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left-width:1px;border-left-color:rgb(204,204,204);border-left-style:solid;padding-left:1ex">
3) What is a local row when I call MatSeqSBAIJSetPreallocationCSR?<br></blockquote><div><br></div><div>Its a row of blocks.</div><div> </div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left-width:1px;border-left-color:rgb(204,204,204);border-left-style:solid;padding-left:1ex">
4) Can I consider my matrix as a single block when calling<br>
MatSeqSBAIJSetPreallocationCSR ?<br></blockquote><div><br></div><div>The maximum block size is 13 I think.</div><div> </div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left-width:1px;border-left-color:rgb(204,204,204);border-left-style:solid;padding-left:1ex">
5) Can I put my acsr array (with values of nonzero elements of the matrix)<br>
as the argument v[] in the MatSeqSBAIJSetPreallocationCSR? The manual<br>
pages told me that v[] should be of the form v[nnz][bs][bs] and I don't<br>
get why not only v[nnz]?<br></blockquote><div><br></div><div>There is no v[] argument:</div><div><br></div><div> <a href="http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatSeqSBAIJSetPreallocation.html">http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatSeqSBAIJSetPreallocation.html</a></div><div><br></div><div> Thanks,</div><div><br></div><div> Matt</div><div> </div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left-width:1px;border-left-color:rgb(204,204,204);border-left-style:solid;padding-left:1ex">
Sorry for that set of questions but I missed an example of usage for a<br>
symmetric CSR matrix.<br>
<br>
Maybe a few lines of code for getting a symmetric (with upper part only<br>
stored) matrix from ia, ja and a arrays (usual CSR format) would help me<br>
to get the whole idea of PETSC matrix creating.<br>
<br>
Thanks in advance!<br>
<br>
<br>
Best regards,<br>
<br>
Kirill Voronin<br>
<br>
<br>
> On Tue, Nov 25, 2014 at 9:07 AM, Kirill Voronin<br>
> <<a href="mailto:kvoronin@labchem.sscc.ru">kvoronin@labchem.sscc.ru</a>><br>
> wrote:<br>
><br>
><br>
>><br>
>> Hello!<br>
>><br>
>><br>
>> I am trying to create a matrix for solving Ax = b with PETSC with given<br>
>> 1-based(!) indexed ia, ja and acsr arrays for a symmetric (only upper<br>
>> part is stored) complex matrix.<br>
>><br>
>> Without memory preallocation it takes too long to fill in all the<br>
>> values of acsr into my Mat A object.<br>
>><br>
>> But I failed to call<br>
>> ierr = MatSeqSBAIJSetPreallocationCSR(Afromcsr, dim, iaHSS, jaHSS,<br>
>> acHSS); because of the 1 based indexing. And I didn't find an option how<br>
>> to handle this issue.<br>
>><br>
>> From the manual:<br>
>> "By default the internal data representation for the AIJ formats employs<br>
>> zero-based indexing. For compatibility with standard Fortran storage,<br>
>> thus enabling use of external Fortran software packages such as SPARSKIT,<br>
>> the option -mat_aij_oneindex enables one-based indexing, where the<br>
>> stored row and column indices begin at one, not zero. All user calls to<br>
>> PETSc routines,<br>
>> regardless of this option, use zero-based indexing."<br>
>><br>
>> So, what can you recommend to do? I obviously don't want to double the<br>
>> memory for storing ia and ja arrays in 0-based indexing. I don't want to<br>
>> change the indexing everywhere to 0-based since I use a piece of<br>
>> external code which works only for 1-based indexing as a preconditioner<br>
>> in PETSC KSP solver.<br>
>><br>
>><br>
><br>
> Just decrement ja, make the call, and increment ja.<br>
><br>
><br>
> Thanks,<br>
><br>
><br>
> Matt<br>
><br>
><br>
><br>
>> Thanks in advance!<br>
>><br>
>><br>
>> With best regards,<br>
>><br>
>><br>
>> Kirill Voronin<br>
>><br>
>><br>
>><br>
>><br>
>><br>
><br>
<span class=""><font color="#888888">><br>
> --<br>
> What most experimenters take for granted before they begin their<br>
> experiments is infinitely more interesting than any results to which their<br>
> experiments lead. -- Norbert Wiener<br>
><br>
><br>
<br>
<br>
--<br>
<br>
<br>
</font></span></blockquote></div><br><br clear="all"><div><br></div>-- <br><div class="gmail_signature">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></div>