[petsc-dev] MatAXPY with DIFFERENT_NONZERO_PATTERN
Alexander Grayver
agrayver at gfz-potsdam.de
Sun Apr 22 08:11:33 CDT 2012
Matt,
Here are matrices (double complex, 3 MB each):
http://dl.dropbox.com/u/60982984/X.dat
http://dl.dropbox.com/u/60982984/Y.dat
Code:
static char help[] = "";
int main(int argc,char **args)
{
Mat A,B;
PetscErrorCode ierr;
PetscViewer viewer;
PetscScalar cone;
PetscInitialize(&argc,&args,(char *)0,help);
ierr =
PetscViewerBinaryOpen(PETSC_COMM_WORLD,"Y.dat",FILE_MODE_READ,&viewer);CHKERRQ(ierr);
ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr);
ierr = MatSetFromOptions(A);CHKERRQ(ierr);
ierr = MatLoad(A,viewer);CHKERRQ(ierr);
ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
ierr =
PetscViewerBinaryOpen(PETSC_COMM_WORLD,"X.dat",FILE_MODE_READ,&viewer);CHKERRQ(ierr);
ierr = MatCreate(PETSC_COMM_WORLD,&B);CHKERRQ(ierr);
ierr = MatSetFromOptions(B);CHKERRQ(ierr);
ierr = MatLoad(B,viewer);CHKERRQ(ierr);
ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
cone = 1.0;
ierr = MatAYPX(A,cone,B,DIFFERENT_NONZERO_PATTERN);
ierr = PetscFinalize();
return 0;
}
It works fine if you run it with:
mpirun -n 10 solveTest -mat_type aij
and crashes if n > 10
Thanks.
On 22.04.2012 14:51, Matthew Knepley wrote:
> On Sun, Apr 22, 2012 at 8:41 AM, Alexander Grayver
> <agrayver at gfz-potsdam.de <mailto:agrayver at gfz-potsdam.de>> wrote:
>
> Sorry, I don't quite get. Is this something I can avoid by doing
> preallocation properly or petsc issue that will be fixed?
>
>
> This looks like a PETSc bug to me, but it would be nice to get input
> that triggers it.
>
> Thanks,
>
> Matt
>
> On 20.04.2012 14:23, Jed Brown wrote:
>> On Fri, Apr 20, 2012 at 05:18, Matthew Knepley <knepley at gmail.com
>> <mailto:knepley at gmail.com>> wrote:
>>
>> Did someone intend to fix this, and then
>> stop? MatAXPY_BasicWithPreallocation() is the same
>> as MatAXPY_Basic().
>>
>>
>> Look how it's called. I don't know what's going on here.
>>
>> ierr =
>> MatAXPYGetPreallocation_SeqAIJ(yy->A,xx->A,nnz_d);CHKERRQ(ierr);
>> ierr =
>> MatAXPYGetPreallocation_SeqAIJ(yy->B,xx->B,nnz_o);CHKERRQ(ierr);
>> ierr =
>> MatMPIAIJSetPreallocation(B,0,nnz_d,0,nnz_o);CHKERRQ(ierr);
>> ierr = MatAXPY_BasicWithPreallocation(B,Y,a,X,str);CHKERRQ(ierr);
>>
>
>
> --
> Regards,
> Alexander
>
>
>
>
> --
> What most experimenters take for granted before they begin their
> experiments is infinitely more interesting than any results to which
> their experiments lead.
> -- Norbert Wiener
--
Regards,
Alexander
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-dev/attachments/20120422/f09cfae0/attachment.html>
More information about the petsc-dev
mailing list