[petsc-users] Matrix Assembly too slow

Vijay Gopal Chilkuri vijay.gopal.c at gmail.com
Sat May 2 17:58:37 CDT 2015


You are right ! I am using the IRPF90 tool to write my program in FORTRAN.

To use IRPF90 you can clone the github repository and get the irpf90 binary
from the bin directory.
The code to generate the library that contains UNIT_L1_ can be found at:
https://github.com/v1j4y/slepc_version.git
just clone this repository and do:

*make*
*make irpf90.a*

That's it !

I'm resending you the PETSc files along with the makefile. you'll need the
irpf90.a library to make the *ex1 *executable.

Hopefully this time it'll work.
Please let me know if something goes wrong.

On Sun, May 3, 2015 at 12:24 AM, Barry Smith <bsmith at mcs.anl.gov> wrote:

>
>   Our stupid mail system blocked the attachment because it contained
> irpf90.a   Is this a library generated with
> https://github.com/scemama/irpf90 ?
> You need to email everything that goes into making that library because it
> is unlikely your .a library will work on my machine. What is this thing and
> can you just tar up the directory that builds it and I can generate irpf90.a
>
>   Barry
>
> > On May 2, 2015, at 3:49 PM, Vijay Gopal Chilkuri <
> vijay.gopal.c at gmail.com> wrote:
> >
> > Ok so here goes.
> >
> > I've attached a tarball with a directory which contains the required
> files to compile ex1.c
> > The UNIT_L1_ subroutine generates the nonzero elements of my hamiltonian
> given a row, it is contained in the
> > irpf90.a library. The actual fortran code is quite ugly.
> > Please tell me if you are able to compile and run the codes.
> >
> > Thanks a lot,
> >  Vijay
> >
> > On Sat, May 2, 2015 at 10:26 PM, Barry Smith <bsmith at mcs.anl.gov> wrote:
> >
> > > On May 2, 2015, at 2:57 PM, Vijay Gopal Chilkuri <
> vijay.gopal.c at gmail.com> wrote:
> > >
> > > Please let me clarify,
> > >
> > > @Matt:
> > > I've followed the proper manner to rapid assembly.
> > > That is :
> > >   1. preallocation
> > >   2. nonzero ownership
> > > Despite of having done the proper preallocation (as you can see from
> the output file),
> > > the matrix takes ages to assemble, whereas the same code works for a
> larger matrix with (100,000,000) elements taking 2 min to assemble !
> > >
> > > @Barry:
> > > I am not touching any file in the code that I attachd.
> > > What i'm doing is that each processor will (according to it's
> ownership range) get it's nonzero entries
> > > using a Fortran subroutine (UNIT_L1_) and call MatSetValues() for
> their respective rows.
> > > Again to be clear I"m NOT doing file I/O anywhere.
> >
> >   Oh sorry, where is the fortran code? I assumed the fortran code was
> reading from a file.
> >
> >    Based on what you sent yes it should take little time; if you send
> something I could run I could try to see why it would be slow.
> >
> >   Barry
> >
> > >
> > > Could you please have a look at the files attached.
> > >
> > >
> > >
> > >
> > >
> > > On Sat, May 2, 2015 at 9:42 PM, Barry Smith <bsmith at mcs.anl.gov>
> wrote:
> > >
> > >   You have two hundred processes banging on the same file on the same
> filesystem fighting with each other to read in the matrix. That's like
> hiring two hundred guys to dig a ditch but making them share one shovel;
> you won't get speed up, you'll get speed down.
> > >
> > >
> > >   Make a simple sequential PETSc program that reads in the matrix and
> then store it with MatView() in binary format. Then in your parallel
> program use MatLoad() to read the matrix in (it will take a few seconds to
> read in the matrix in that case).
> > >
> > >   Barry
> > >
> > >
> > > > On May 2, 2015, at 2:10 PM, Vijay Gopal Chilkuri <
> vijay.gopal.c at gmail.com> wrote:
> > > >
> > > > Hi,
> > > >
> > > > I'm trying to diagonalize large matrices using PETSc and SLEPc.
> > > > I've successfully diagonalized 100 million dimensional matrix in
> 2hours.
> > > >
> > > > But, the problem is when i try to diagonailze (similar) a smaller
> metrix of
> > > > dimension 67,000,000 (24 nonzero elements per row) with 10 noes
> running 20 processors each (similar to the 100 million case), the matrix
> assembly itself
> > > > takes 2h !
> > > >
> > > > Could someone point out the mistakes that i'm making ?
> > > >
> > > > i attach the source code and the output with the mail.
> > > >
> > > > thanks,
> > > >  Vijay
> > > >
> > > >
> > > > #include <slepceps.h>
> > > > #include <petsctime.h>
> > > >
> > > > #undef __FUNCT__
> > > > #define __FUNCT__ "main"
> > > > int main(int argc,char **argv)
> > > > {
> > > >   Mat            A;           /* problem matrix */
> > > >   EPS            eps;         /* eigenproblem solver context */
> > > >   EPSType        type;
> > > >   PetscReal      error,tol,re,im;
> > > >   PetscScalar    kr,ki,value[52];
> > > >   Vec            xr,xi;
> > > >   PetscInt
>  n=16224936,ii,i,veclen,j,Istart,Iend,col[52],nev,maxit,its,nconv,countcol;
> > > >   PetscInt       d_nz,o_nz;
> > > >   PetscLogDouble t1,t2,tt1,tt2;
> > > >   PetscBool      FirstBlock=PETSC_FALSE,LastBlock=PETSC_FALSE;
> > > >   PetscErrorCode ierr;
> > > >   PetscScalar    eigr;
> > > >   PetscScalar    eigi;
> > > >   PetscScalar    * data;
> > > >   Vec          Vr,Vi;
> > > >   char           filename[PETSC_MAX_PATH_LEN]="FIL666";
> > > >   PetscViewer    viewer;
> > > >   PetscBool      ishermitian;
> > > >   int            mpiid;
> > > >   long int            kk,iii;
> > > >   long int            tcountcol,tcol[52];
> > > >   float          val[52];
> > > >   long int            ntrou=1;
> > > >   long int            l1[52]={1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 1,
> 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17, 18, 19, 20, 21, 22,
> 23, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0       };
> > > >   long int            l2[52]={2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12,
> 24, 23, 22, 21, 20, 19, 18, 17, 16, 15, 14, 13, 14, 15, 16, 17, 18, 19, 20,
> 21, 22, 23, 24, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0       };
> > > >   long int          ktyp[52]={1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 2, 2,
> 2, 2, 2, 2, 2, 2, 2, 2, 2, 2, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 3, 0, 0, 0, 0,
> 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0      };
> > > >   double           xjjz[52]  ={0.0333333333333,-0.8, 0.,
> > > >       0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0,
> 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0};
> > > >   double           xjjxy[52] ={0.0333333333333,-0.8,0.,
> > > >       0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0,
> 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0};
> > > >   double           xtt[52]   ={-1.0,0.,0.,
> > > >       0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0,
> 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0};
> > > >
> > > >   SlepcInitialize(&argc,&argv,(char*)0,NULL);
> > > >   ierr = PetscPrintf(PETSC_COMM_WORLD,"\n1-D t-J Eigenproblem,
> n=%D\n\n",n);CHKERRQ(ierr);
> > > >   ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr);
> > > > //ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);CHKERRQ(ierr);
> > > >   ierr =
> MatCreateAIJ(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,n,n,46,NULL,46,NULL,&A);CHKERRQ(ierr);
> > > >   ierr = MatSetFromOptions(A);CHKERRQ(ierr);
> > > >   ierr = MatSetUp(A);CHKERRQ(ierr);
> > > >
> > > >   ierr = MatGetOwnershipRange(A,&Istart,&Iend);CHKERRQ(ierr);
> > > >   ierr =
> MatCreateAIJ(PETSC_COMM_WORLD,Iend-Istart,Iend-Istart,n,n,46,NULL,46,NULL,&A);CHKERRQ(ierr);
> > > >   MPI_Comm_rank(MPI_COMM_WORLD,&mpiid);
> > > >   ierr = PetscTime(&tt1);CHKERRQ(ierr);
> > > >   for (i=Istart; i<Iend; i++) {
> > > >       tcountcol=0;
> > > >       iii=i+1;
> > > >       if(i%5 == 0 && mpiid==0){
> > > >         ierr = PetscTime(&t1);CHKERRQ(ierr);
> > > >       }
> > > >       unit_l1_(
> > > >             l1,
> > > >             l2,
> > > >             ktyp,
> > > >             &iii,
> > > >             xjjxy,
> > > >             xjjz ,
> > > >             xtt  ,
> > > >             &tcountcol,
> > > >             &ntrou,
> > > >             tcol,
> > > >             val);
> > > >       for(kk=0;kk<52;kk++){
> > > >           value[kk] = val[kk];
> > > >           col[kk] = tcol[kk]-1;
> > > >       }
> > > >       if(i%5 == 0 && mpiid==0){
> > > >         ierr = PetscTime(&t2);CHKERRQ(ierr);
> > > >         ierr = PetscPrintf(PETSC_COMM_WORLD," i: %d\n mpiid:
> %d\ntime: %f\n",i,mpiid,t2-t1);CHKERRQ(ierr);
> > > >       }
> > > >       countcol=tcountcol;
> > > >       if(i%5 == 0 && mpiid==0){
> > > >         ierr = PetscTime(&t1);CHKERRQ(ierr);
> > > >       }
> > > >     ierr =
> MatSetValues(A,1,&i,countcol,col,value,INSERT_VALUES);CHKERRQ(ierr);
> > > >       if(i%5 == 0 && mpiid==0){
> > > >         ierr = PetscTime(&t2);CHKERRQ(ierr);
> > > >         ierr = PetscPrintf(PETSC_COMM_WORLD," processor \n mpiid:
> %d\ntime: %f\n",mpiid,t2-t1);CHKERRQ(ierr);
> > > >       }
> > > >   }
> > > >   ierr = PetscTime(&tt2);CHKERRQ(ierr);
> > > >   ierr = PetscPrintf(PETSC_COMM_WORLD," Time used to build the
> matrix: %f\n",tt2-tt1);CHKERRQ(ierr);
> > > >
> > > >
> > > >   ierr = PetscTime(&tt1);CHKERRQ(ierr);
> > > >   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
> > > >   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
> > > >   ierr = PetscTime(&tt2);CHKERRQ(ierr);
> > > >   ierr = PetscPrintf(PETSC_COMM_WORLD," Time used to assemble the
> matrix: %f\n",tt2-tt1);CHKERRQ(ierr);
> > > >   ierr = MatGetVecs(A,NULL,&xr);CHKERRQ(ierr);
> > > >   ierr = MatGetVecs(A,NULL,&xi);CHKERRQ(ierr);
> > > >
> > > >   ierr = EPSCreate(PETSC_COMM_WORLD,&eps);CHKERRQ(ierr);
> > > >   ierr = EPSSetOperators(eps,A,NULL);CHKERRQ(ierr);
> > > >   ierr = EPSSetProblemType(eps,EPS_HEP);CHKERRQ(ierr);
> > > >   ierr = EPSSetWhichEigenpairs(eps,EPS_SMALLEST_REAL);CHKERRQ(ierr);
> > > >
> > > >   ierr = EPSSetFromOptions(eps);CHKERRQ(ierr);
> > > >   tol = 1.e-8;
> > > >   maxit = 10000000;
> > > >   ierr = EPSSetTolerances(eps,tol,maxit);CHKERRQ(ierr);
> > > >
> > > >   ierr = PetscTime(&t1);CHKERRQ(ierr);
> > > >   ierr = EPSSolve(eps);CHKERRQ(ierr);
> > > >   ierr = PetscTime(&t2);CHKERRQ(ierr);
> > > >   ierr = PetscPrintf(PETSC_COMM_WORLD," Time used:
> %f\n",t2-t1);CHKERRQ(ierr);
> > > >   ierr = EPSGetIterationNumber(eps,&its);CHKERRQ(ierr);
> > > >   ierr = PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the
> method: %D\n",its);CHKERRQ(ierr);
> > > >   ierr = EPSGetType(eps,&type);CHKERRQ(ierr);
> > > >   ierr = PetscPrintf(PETSC_COMM_WORLD," Solution method:
> %s\n\n",type);CHKERRQ(ierr);
> > > >   ierr = EPSGetDimensions(eps,&nev,NULL,NULL);CHKERRQ(ierr);
> > > >   ierr = PetscPrintf(PETSC_COMM_WORLD," Number of requested
> eigenvalues: %D\n",nev);CHKERRQ(ierr);
> > > >   ierr = EPSGetTolerances(eps,&tol,&maxit);CHKERRQ(ierr);
> > > >   ierr = PetscPrintf(PETSC_COMM_WORLD," Stopping condition:
> tol=%.4g, maxit=%D\n",(double)tol,maxit);CHKERRQ(ierr);
> > > >
> > > >   ierr = EPSGetConverged(eps,&nconv);CHKERRQ(ierr);
> > > >   ierr = EPSPrintSolution(eps,NULL);CHKERRQ(ierr);
> > > >   /*
> > > >   EPSGetConverged(eps,&nconv);
> > > >   if (nconv>0) {
> > > >     PetscViewerASCIIOpen(PETSC_COMM_WORLD,filename,&viewer);
> > > >     EPSIsHermitian(eps,&ishermitian);
> > > >     for (i=0;i<nconv;i++) {
> > > >       EPSGetEigenvector(eps,i,xr,xi);
> > > >       VecView(xr,viewer);
> > > > #if !defined(PETSC_USE_COMPLEX)
> > > >       if (!ishermitian) { VecView(xi,viewer); }
> > > > #endif
> > > >     }
> > > >     PetscViewerDestroy(&viewer);
> > > >   }
> > > >   */
> > > >   ierr = EPSDestroy(&eps);CHKERRQ(ierr);
> > > >   ierr = MatDestroy(&A);CHKERRQ(ierr);
> > > >   ierr = VecDestroy(&xr);CHKERRQ(ierr);
> > > >   ierr = VecDestroy(&xi);CHKERRQ(ierr);
> > > >   ierr = SlepcFinalize();
> > > >   return 0;
> > > > }
> > > >
> > > > OUTPUT:
> > > >  Time used to build the matrix: 1914.729022
> > > > [0] MatStashScatterBegin_Private(): No of messages: 0
> > > > [0] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [6] MatStashScatterBegin_Private(): No of messages: 0
> > > > [6] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [182] MatStashScatterBegin_Private(): No of messages: 0
> > > > [182] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [100] MatStashScatterBegin_Private(): No of messages: 0
> > > > [100] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [60] MatStashScatterBegin_Private(): No of messages: 0
> > > > [60] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [140] MatStashScatterBegin_Private(): No of messages: 0
> > > > [140] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [120] MatStashScatterBegin_Private(): No of messages: 0
> > > > [120] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [162] MatStashScatterBegin_Private(): No of messages: 0
> > > > [162] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [20] MatStashScatterBegin_Private(): No of messages: 0
> > > > [20] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [80] MatStashScatterBegin_Private(): No of messages: 0
> > > > [80] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [40] MatStashScatterBegin_Private(): No of messages: 0
> > > > [40] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [7] MatStashScatterBegin_Private(): No of messages: 0
> > > > [7] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [185] MatStashScatterBegin_Private(): No of messages: 0
> > > > [185] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [101] MatStashScatterBegin_Private(): No of messages: 0
> > > > [101] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [61] MatStashScatterBegin_Private(): No of messages: 0
> > > > [61] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [142] MatStashScatterBegin_Private(): No of messages: 0
> > > > [142] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [130] MatStashScatterBegin_Private(): No of messages: 0
> > > > [130] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [166] MatStashScatterBegin_Private(): No of messages: 0
> > > > [166] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [21] MatStashScatterBegin_Private(): No of messages: 0
> > > > [21] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [81] MatStashScatterBegin_Private(): No of messages: 0
> > > > [81] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [44] MatStashScatterBegin_Private(): No of messages: 0
> > > > [44] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [9] MatStashScatterBegin_Private(): No of messages: 0
> > > > [9] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [187] MatStashScatterBegin_Private(): No of messages: 0
> > > > [187] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [105] MatStashScatterBegin_Private(): No of messages: 0
> > > > [105] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [62] MatStashScatterBegin_Private(): No of messages: 0
> > > > [62] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [143] MatStashScatterBegin_Private(): No of messages: 0
> > > > [143] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [132] MatStashScatterBegin_Private(): No of messages: 0
> > > > [132] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [169] MatStashScatterBegin_Private(): No of messages: 0
> > > > [169] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [24] MatStashScatterBegin_Private(): No of messages: 0
> > > > [24] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [84] MatStashScatterBegin_Private(): No of messages: 0
> > > > [84] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [45] MatStashScatterBegin_Private(): No of messages: 0
> > > > [45] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [10] MatStashScatterBegin_Private(): No of messages: 0
> > > > [10] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [188] MatStashScatterBegin_Private(): No of messages: 0
> > > > [188] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [109] MatStashScatterBegin_Private(): No of messages: 0
> > > > [109] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [63] MatStashScatterBegin_Private(): No of messages: 0
> > > > [63] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [148] MatStashScatterBegin_Private(): No of messages: 0
> > > > [148] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [133] MatStashScatterBegin_Private(): No of messages: 0
> > > > [133] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [170] MatStashScatterBegin_Private(): No of messages: 0
> > > > [170] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [25] MatStashScatterBegin_Private(): No of messages: 0
> > > > [25] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [86] MatStashScatterBegin_Private(): No of messages: 0
> > > > [86] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [47] MatStashScatterBegin_Private(): No of messages: 0
> > > > [47] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [11] MatStashScatterBegin_Private(): No of messages: 0
> > > > [11] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [190] MatStashScatterBegin_Private(): No of messages: 0
> > > > [190] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [111] MatStashScatterBegin_Private(): No of messages: 0
> > > > [111] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [64] MatStashScatterBegin_Private(): No of messages: 0
> > > > [64] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [156] MatStashScatterBegin_Private(): No of messages: 0
> > > > [156] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [134] MatStashScatterBegin_Private(): No of messages: 0
> > > > [134] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [171] MatStashScatterBegin_Private(): No of messages: 0
> > > > [171] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [26] MatStashScatterBegin_Private(): No of messages: 0
> > > > [26] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [91] MatStashScatterBegin_Private(): No of messages: 0
> > > > [91] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [48] MatStashScatterBegin_Private(): No of messages: 0
> > > > [48] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [12] MatStashScatterBegin_Private(): No of messages: 0
> > > > [12] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [192] MatStashScatterBegin_Private(): No of messages: 0
> > > > [192] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [112] MatStashScatterBegin_Private(): No of messages: 0
> > > > [112] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [65] MatStashScatterBegin_Private(): No of messages: 0
> > > > [65] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [157] MatStashScatterBegin_Private(): No of messages: 0
> > > > [157] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [135] MatStashScatterBegin_Private(): No of messages: 0
> > > > [135] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [174] MatStashScatterBegin_Private(): No of messages: 0
> > > > [174] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [28] MatStashScatterBegin_Private(): No of messages: 0
> > > > [28] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [13] MatStashScatterBegin_Private(): No of messages: 0
> > > > [13] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [194] MatStashScatterBegin_Private(): No of messages: 0
> > > > [194] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [113] MatStashScatterBegin_Private(): No of messages: 0
> > > > [113] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [66] MatStashScatterBegin_Private(): No of messages: 0
> > > > [66] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [158] MatStashScatterBegin_Private(): No of messages: 0
> > > > [158] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [136] MatStashScatterBegin_Private(): No of messages: 0
> > > > [136] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [160] MatStashScatterBegin_Private(): No of messages: 0
> > > > [160] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [29] MatStashScatterBegin_Private(): No of messages: 0
> > > > [29] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [92] MatStashScatterBegin_Private(): No of messages: 0
> > > > [92] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [53] MatStashScatterBegin_Private(): No of messages: 0
> > > > [53] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [14] MatStashScatterBegin_Private(): No of messages: 0
> > > > [14] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [195] MatStashScatterBegin_Private(): No of messages: 0
> > > > [195] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [114] MatStashScatterBegin_Private(): No of messages: 0
> > > > [114] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [67] MatStashScatterBegin_Private(): No of messages: 0
> > > > [67] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [141] MatStashScatterBegin_Private(): No of messages: 0
> > > > [141] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [137] MatStashScatterBegin_Private(): No of messages: 0
> > > > [137] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [161] MatStashScatterBegin_Private(): No of messages: 0
> > > > [161] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [32] MatStashScatterBegin_Private(): No of messages: 0
> > > > [32] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [93] MatStashScatterBegin_Private(): No of messages: 0
> > > > [93] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [54] MatStashScatterBegin_Private(): No of messages: 0
> > > > [54] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [15] MatStashScatterBegin_Private(): No of messages: 0
> > > > [15] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [196] MatStashScatterBegin_Private(): No of messages: 0
> > > > [196] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [115] MatStashScatterBegin_Private(): No of messages: 0
> > > > [115] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [68] MatStashScatterBegin_Private(): No of messages: 0
> > > > [68] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [144] MatStashScatterBegin_Private(): No of messages: 0
> > > > [144] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [138] MatStashScatterBegin_Private(): No of messages: 0
> > > > [138] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [163] MatStashScatterBegin_Private(): No of messages: 0
> > > > [163] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [33] MatStashScatterBegin_Private(): No of messages: 0
> > > > [33] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [94] MatStashScatterBegin_Private(): No of messages: 0
> > > > [94] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [55] MatStashScatterBegin_Private(): No of messages: 0
> > > > [55] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [17] MatStashScatterBegin_Private(): No of messages: 0
> > > > [17] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [198] MatStashScatterBegin_Private(): No of messages: 0
> > > > [198] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [116] MatStashScatterBegin_Private(): No of messages: 0
> > > > [116] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [69] MatStashScatterBegin_Private(): No of messages: 0
> > > > [69] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [145] MatStashScatterBegin_Private(): No of messages: 0
> > > > [145] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [139] MatStashScatterBegin_Private(): No of messages: 0
> > > > [139] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [164] MatStashScatterBegin_Private(): No of messages: 0
> > > > [164] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [34] MatStashScatterBegin_Private(): No of messages: 0
> > > > [34] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [95] MatStashScatterBegin_Private(): No of messages: 0
> > > > [95] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [56] MatStashScatterBegin_Private(): No of messages: 0
> > > > [56] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [1] MatStashScatterBegin_Private(): No of messages: 0
> > > > [1] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [199] MatStashScatterBegin_Private(): No of messages: 0
> > > > [199] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [119] MatStashScatterBegin_Private(): No of messages: 0
> > > > [119] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [70] MatStashScatterBegin_Private(): No of messages: 0
> > > > [70] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [146] MatStashScatterBegin_Private(): No of messages: 0
> > > > [146] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [121] MatStashScatterBegin_Private(): No of messages: 0
> > > > [121] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [165] MatStashScatterBegin_Private(): No of messages: 0
> > > > [165] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [35] MatStashScatterBegin_Private(): No of messages: 0
> > > > [35] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [99] MatStashScatterBegin_Private(): No of messages: 0
> > > > [99] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [59] MatStashScatterBegin_Private(): No of messages: 0
> > > > [59] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [2] MatStashScatterBegin_Private(): No of messages: 0
> > > > [2] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [180] MatStashScatterBegin_Private(): No of messages: 0
> > > > [180] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [102] MatStashScatterBegin_Private(): No of messages: 0
> > > > [102] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [71] MatStashScatterBegin_Private(): No of messages: 0
> > > > [71] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [147] MatStashScatterBegin_Private(): No of messages: 0
> > > > [147] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [122] MatStashScatterBegin_Private(): No of messages: 0
> > > > [122] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [167] MatStashScatterBegin_Private(): No of messages: 0
> > > > [167] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [36] MatStashScatterBegin_Private(): No of messages: 0
> > > > [36] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [82] MatStashScatterBegin_Private(): No of messages: 0
> > > > [82] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [41] MatStashScatterBegin_Private(): No of messages: 0
> > > > [41] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [3] MatStashScatterBegin_Private(): No of messages: 0
> > > > [3] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [181] MatStashScatterBegin_Private(): No of messages: 0
> > > > [181] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [103] MatStashScatterBegin_Private(): No of messages: 0
> > > > [103] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [72] MatStashScatterBegin_Private(): No of messages: 0
> > > > [72] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [149] MatStashScatterBegin_Private(): No of messages: 0
> > > > [149] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [123] MatStashScatterBegin_Private(): No of messages: 0
> > > > [123] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [168] MatStashScatterBegin_Private(): No of messages: 0
> > > > [168] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [37] MatStashScatterBegin_Private(): No of messages: 0
> > > > [37] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [83] MatStashScatterBegin_Private(): No of messages: 0
> > > > [83] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [42] MatStashScatterBegin_Private(): No of messages: 0
> > > > [42] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [4] MatStashScatterBegin_Private(): No of messages: 0
> > > > [4] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [183] MatStashScatterBegin_Private(): No of messages: 0
> > > > [183] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [104] MatStashScatterBegin_Private(): No of messages: 0
> > > > [104] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [73] MatStashScatterBegin_Private(): No of messages: 0
> > > > [73] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [150] MatStashScatterBegin_Private(): No of messages: 0
> > > > [150] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [124] MatStashScatterBegin_Private(): No of messages: 0
> > > > [124] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [172] MatStashScatterBegin_Private(): No of messages: 0
> > > > [172] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [38] MatStashScatterBegin_Private(): No of messages: 0
> > > > [38] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [85] MatStashScatterBegin_Private(): No of messages: 0
> > > > [85] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [43] MatStashScatterBegin_Private(): No of messages: 0
> > > > [43] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [5] MatStashScatterBegin_Private(): No of messages: 0
> > > > [5] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [184] MatStashScatterBegin_Private(): No of messages: 0
> > > > [184] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [106] MatStashScatterBegin_Private(): No of messages: 0
> > > > [106] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [74] MatStashScatterBegin_Private(): No of messages: 0
> > > > [74] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [151] MatStashScatterBegin_Private(): No of messages: 0
> > > > [151] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [125] MatStashScatterBegin_Private(): No of messages: 0
> > > > [125] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [173] MatStashScatterBegin_Private(): No of messages: 0
> > > > [173] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [39] MatStashScatterBegin_Private(): No of messages: 0
> > > > [39] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [87] MatStashScatterBegin_Private(): No of messages: 0
> > > > [87] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [46] MatStashScatterBegin_Private(): No of messages: 0
> > > > [46] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [8] MatStashScatterBegin_Private(): No of messages: 0
> > > > [8] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [186] MatStashScatterBegin_Private(): No of messages: 0
> > > > [186] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [107] MatStashScatterBegin_Private(): No of messages: 0
> > > > [107] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [75] MatStashScatterBegin_Private(): No of messages: 0
> > > > [75] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [152] MatStashScatterBegin_Private(): No of messages: 0
> > > > [152] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [126] MatStashScatterBegin_Private(): No of messages: 0
> > > > [126] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [175] MatStashScatterBegin_Private(): No of messages: 0
> > > > [175] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [22] MatStashScatterBegin_Private(): No of messages: 0
> > > > [22] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [88] MatStashScatterBegin_Private(): No of messages: 0
> > > > [88] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [49] MatStashScatterBegin_Private(): No of messages: 0
> > > > [49] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [16] MatStashScatterBegin_Private(): No of messages: 0
> > > > [16] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [189] MatStashScatterBegin_Private(): No of messages: 0
> > > > [189] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [108] MatStashScatterBegin_Private(): No of messages: 0
> > > > [108] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [76] MatStashScatterBegin_Private(): No of messages: 0
> > > > [76] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [153] MatStashScatterBegin_Private(): No of messages: 0
> > > > [153] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [127] MatStashScatterBegin_Private(): No of messages: 0
> > > > [127] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [176] MatStashScatterBegin_Private(): No of messages: 0
> > > > [176] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [23] MatStashScatterBegin_Private(): No of messages: 0
> > > > [23] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [89] MatStashScatterBegin_Private(): No of messages: 0
> > > > [89] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [50] MatStashScatterBegin_Private(): No of messages: 0
> > > > [50] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [18] MatStashScatterBegin_Private(): No of messages: 0
> > > > [18] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [191] MatStashScatterBegin_Private(): No of messages: 0
> > > > [191] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [110] MatStashScatterBegin_Private(): No of messages: 0
> > > > [110] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [77] MatStashScatterBegin_Private(): No of messages: 0
> > > > [77] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [154] MatStashScatterBegin_Private(): No of messages: 0
> > > > [154] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [128] MatStashScatterBegin_Private(): No of messages: 0
> > > > [128] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [177] MatStashScatterBegin_Private(): No of messages: 0
> > > > [177] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [27] MatStashScatterBegin_Private(): No of messages: 0
> > > > [27] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [90] MatStashScatterBegin_Private(): No of messages: 0
> > > > [90] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [51] MatStashScatterBegin_Private(): No of messages: 0
> > > > [51] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [19] MatStashScatterBegin_Private(): No of messages: 0
> > > > [19] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [193] MatStashScatterBegin_Private(): No of messages: 0
> > > > [193] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [117] MatStashScatterBegin_Private(): No of messages: 0
> > > > [117] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [78] MatStashScatterBegin_Private(): No of messages: 0
> > > > [78] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [155] MatStashScatterBegin_Private(): No of messages: 0
> > > > [155] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [129] MatStashScatterBegin_Private(): No of messages: 0
> > > > [129] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [178] MatStashScatterBegin_Private(): No of messages: 0
> > > > [178] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [30] MatStashScatterBegin_Private(): No of messages: 0
> > > > [30] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [96] MatStashScatterBegin_Private(): No of messages: 0
> > > > [96] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [52] MatStashScatterBegin_Private(): No of messages: 0
> > > > [52] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [197] MatStashScatterBegin_Private(): No of messages: 0
> > > > [197] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [118] MatStashScatterBegin_Private(): No of messages: 0
> > > > [118] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [79] MatStashScatterBegin_Private(): No of messages: 0
> > > > [79] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [159] MatStashScatterBegin_Private(): No of messages: 0
> > > > [159] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [131] MatStashScatterBegin_Private(): No of messages: 0
> > > > [131] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [179] MatStashScatterBegin_Private(): No of messages: 0
> > > > [179] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [31] MatStashScatterBegin_Private(): No of messages: 0
> > > > [31] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [97] MatStashScatterBegin_Private(): No of messages: 0
> > > > [97] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [57] MatStashScatterBegin_Private(): No of messages: 0
> > > > [57] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [98] MatStashScatterBegin_Private(): No of messages: 0
> > > > [98] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [58] MatStashScatterBegin_Private(): No of messages: 0
> > > > [58] MatAssemblyBegin_MPIAIJ(): Stash has 0 entries, uses 0 mallocs.
> > > > [104] MatAssemblyEnd_SeqAIJ(): Matrix size: 81125 X 81125; storage
> space: 3238253 unneeded,493497 used
> > > > [104] MatAssemblyEnd_SeqAIJ(): Number of mallocs during
> MatSetValues() is 0
> > > > [104] MatAssemblyEnd_SeqAIJ(): Maximum nonzeros in any row is 12
> > > > [104] MatCheckCompressedRow(): Found the ratio (num_zerorows
> 0)/(num_localrows 81125) < 0.6. Do not use CompressedRow routines.
> > > > [4] MatAssemblyEnd_SeqAIJ(): Matrix size: 81125 X 81125; storage
> space: 3197999 unneeded,533751 used
> > > > [4] MatAssemblyEnd_SeqAIJ(): Number of mallocs during MatSetValues()
> is 0
> > > > [4] MatAssemblyEnd_SeqAIJ(): Maximum nonzeros in any row is 13
> > > > [4] MatCheckCompressedRow(): Found the ratio (num_zerorows
> 0)/(num_localrows 81125) < 0.6. Do not use CompressedRow routines.
> > > > [106] MatAssemblyEnd_SeqAIJ(): Matrix size: 81125 X 81125; storage
> space: 3243809 unneeded,487941 used
> > > > [106] MatAssemblyEnd_SeqAIJ(): Number of mallocs during
> MatSetValues() is 0
> > > > [106] MatAssemblyEnd_SeqAIJ(): Maximum nonzeros in any row is 12
> > > > [106] MatCheckCompressedRow(): Found the ratio (num_zerorows
> 0)/(num_localrows 81125) < 0.6. Do not use CompressedRow routines.
> > > > [144] MatAssemblyEnd_SeqAIJ(): Matrix size: 81124 X 81124; storage
> space: 3243446 unneeded,488258 used
> > > > [144] MatAssemblyEnd_SeqAIJ(): Number of mallocs during
> MatSetValues() is 0
> > > > [144] MatAssemblyEnd_SeqAIJ(): Maximum nonzeros in any row is 11
> > > > [144] MatCheckCompressedRow(): Found the ratio (num_zerorows
> 0)/(num_localrows 81124) < 0.6. Do not use CompressedRow routines.
> > > > [123] MatAssemblyEnd_SeqAIJ(): Matrix size: 81125 X 81125; storage
> space: 3237359 unneeded,494391 used
> > > > [123] MatAssemblyEnd_SeqAIJ(): Number of mallocs during
> MatSetValues() is 0
> > > > [123] MatAssemblyEnd_SeqAIJ(): Maximum nonzeros in any row is 12
> > > > [123] MatCheckCompressedRow(): Found the ratio (num_zerorows
> 0)/(num_localrows 81125) < 0.6. Do not use CompressedRow routines.
> > > > [173] MatAssemblyEnd_SeqAIJ(): Matrix size: 81124 X 81124; storage
> space: 3245082 unneeded,486622 used
> > > > [173] MatAssemblyEnd_SeqAIJ(): Number of mallocs during
> MatSetValues() is 0
> > > > [173] MatAssemblyEnd_SeqAIJ(): Maximum nonzeros in any row is 12
> > > > [173] MatCheckCompressedRow(): Found the ratio (num_zerorows
> 0)/(num_localrows 81124) < 0.6. Do not use CompressedRow routines.
> > > > [85] MatAssemblyEnd_SeqAIJ(): Matrix size: 81125 X 81125; storage
> space: 3217795 unneeded,513955 used
> > > > [85] MatAssemblyEnd_SeqAIJ(): Number of mallocs during
> MatSetValues() is 0
> > > > [85] MatAssemblyEnd_SeqAIJ(): Maximum nonzeros in any row is 12
> > > > [85] MatCheckCompressedRow(): Found the ratio (num_zerorows
> 0)/(num_localrows 81125) < 0.6. Do not use CompressedRow routines.
> > > > [6] MatAssemblyEnd_SeqAIJ(): Matrix size: 81125 X 81125; storage
> space: 3208633 unneeded,523117 used
> > > > [6] MatAssemblyEnd_SeqAIJ(): Number of mallocs during MatSetValues()
> is 0
> > > > [6] MatAssemblyEnd_SeqAIJ(): Maximum nonzeros in any row is 13
> > > > [6] MatCheckCompressedRow(): Found the ratio (num_zerorows
> 0)/(num_localrows 81125) < 0.6. Do not use CompressedRow routines.
> > > > [185] MatAssemblyEnd_SeqAIJ(): Matrix size: 81124 X 81124; storage
> space: 3209020 unneeded,522684 used
> > > > [185] MatAssemblyEnd_SeqAIJ(): Number of mallocs during
> MatSetValues() is 0
> > > > [185] MatAssemblyEnd_SeqAIJ(): Maximum nonzeros in any row is 12
> > > > [185] MatCheckCompressedRow(): Found the ratio (num_zerorows
> 0)/(num_localrows 81124) < 0.6. Do not use CompressedRow routines.
> > > > [108] MatAssemblyEnd_SeqAIJ(): Matrix size: 81125 X 81125; storage
> space: 3216305 unneeded,515445 used
> > > > [108] MatAssemblyEnd_SeqAIJ(): Number of mallocs during
> MatSetValues() is 0
> > > > [108] MatAssemblyEnd_SeqAIJ(): Maximum nonzeros in any row is 12
> > > > [73] MatAssemblyEnd_SeqAIJ(): Matrix size: 81125 X 81125; storage
> space: 3201643 unneeded,530107 used
> > > > [73] MatAssemblyEnd_SeqAIJ(): Number of mallocs during
> MatSetValues() is 0
> > > > [73] MatAssemblyEnd_SeqAIJ(): Maximum nonzeros in any row is 13
> > > > [73] MatCheckCompressedRow(): Found the ratio (num_zerorows
> 0)/(num_localrows 81125) < 0.6. Do not use CompressedRow routines.
> > > > [146] MatAssemblyEnd_SeqAIJ(): Matrix size: 81124 X 81124; storage
> space: 3234040 unneeded,497664 used
> > > > [146] MatAssemblyEnd_SeqAIJ(): Number of mallocs during
> MatSetValues() is 0
> > > > [146] MatAssemblyEnd_SeqAIJ(): Maximum nonzeros in any row is 11
> > > > [146] MatCheckCompressedRow(): Found the ratio (num_zerorows
> 0)/(num_localrows 81124) < 0.6. Do not use CompressedRow routines.
> > > > [127] MatAssemblyEnd_SeqAIJ(): Matrix size: 81125 X 81125; storage
> space: 3243435 unneeded,488315 used
> > > > [127] MatAssemblyEnd_SeqAIJ(): Number of mallocs during
> MatSetValues() is 0
> > > > [127] MatAssemblyEnd_SeqAIJ(): Maximum nonzeros in any row is 12
> > > > [127] MatCheckCompressedRow(): Found the ratio (num_zerorows
> 0)/(num_localrows 81125) < 0.6. Do not use CompressedRow routines.
> > > > [177] MatAssemblyEnd_SeqAIJ(): Matrix size: 81124 X 81124; storage
> space: 3249858 unneeded,481846 used
> > > > [177] MatAssemblyEnd_SeqAIJ(): Number of mallocs during
> MatSetValues() is 0
> > > > [177] MatAssemblyEnd_SeqAIJ(): Maximum nonzeros in any row is 12
> > > > [177] MatCheckCompressedRow(): Found the ratio (num_zerorows
> 0)/(num_localrows 81124) < 0.6. Do not use CompressedRow routines.
> > > > [27] MatAssemblyEnd_SeqAIJ(): Matrix size: 81125 X 81125; storage
> space: 3207727 unneeded,524023 used
> > > > [27] MatAssemblyEnd_SeqAIJ(): Number of mallocs during
> MatSetValues() is 0
> > > > [27] MatAssemblyEnd_SeqAIJ(): Maximum nonzeros in any row is 13
> > > > [27] MatCheckCompressedRow(): Found the ratio (num_zerorows
> 0)/(num_localrows 81125) < 0.6. Do not use CompressedRow routines.
> > > > [87] MatAssemblyEnd_SeqAIJ(): Matrix size: 81125 X 81125; storage
> space: 3220647 unneeded,511103 used
> > > > [87] MatAssemblyEnd_SeqAIJ(): Number of mallocs during
> MatSetValues() is 0
> > > > [87] MatAssemblyEnd_SeqAIJ(): Maximum nonzeros in any row is 12
> > > > [87] MatCheckCompressedRow(): Found the ratio (num_zerorows
> 0)/(num_localrows 81125) < 0.6. Do not use CompressedRow routines.
> > > > [54] MatAssemblyEnd_SeqAIJ(): Matrix size: 81125 X 81125; storage
> space: 3197939 unneeded,533811 used
> > > > [54] MatAssemblyEnd_SeqAIJ(): Number of mallocs during
> MatSetValues() is 0
> > > > [54] MatAssemblyEnd_SeqAIJ(): Maximum nonzeros in any row is 13
> > > > [54] MatCheckCompressedRow(): Found the ratio (num_zerorows
> 0)/(num_localrows 81125) < 0.6. Do not use CompressedRow routines.
> > > > [8] MatAssemblyEnd_SeqAIJ(): Matrix size: 81125 X 81125; storage
> space: 3179441 unneeded,552309 used
> > > > [8] MatAssemblyEnd_SeqAIJ(): Number of mallocs during MatSetValues()
> is 0
> > > > [8] MatAssemblyEnd_SeqAIJ(): Maximum nonzeros in any row is 13
> > > > [8] MatCheckCompressedRow(): Found the ratio (num_zerorows
> 0)/(num_localrows 81125) < 0.6. Do not use CompressedRow routines.
> > > > [187] MatAssemblyEnd_SeqAIJ(): Matrix size: 81124 X 81124; storage
> space: 3197890 unneeded,533814 used
> > > > [187] MatAssemblyEnd_SeqAIJ(): Number of mallocs during
> MatSetValues() is 0
> > > > [187] MatAssemblyEnd_SeqAIJ(): Maximum nonzeros in any row is 13
> > > > [187] MatCheckCompressedRow(): Found the ratio (num_zerorows
> 0)/(num_localrows 81124) < 0.6. Do not use CompressedRow routines.
> > > > [108] MatCheckCompressedRow(): Found the ratio (num_zerorows
> 0)/(num_localrows 81125) < 0.6. Do not use CompressedRow routines.
> > > > [74] MatAssemblyEnd_SeqAIJ(): Matrix size: 81125 X 81125; storage
> space: 3200381 unneeded,531369 used
> > > >
> > > > <ex1.c><inp.log>
> > >
> > >
> >
> >
> > <1_Warning.txt>
>
>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20150503/f133b272/attachment-0001.html>
-------------- next part --------------
A non-text attachment was scrubbed...
Name: ex1.c
Type: text/x-csrc
Size: 6249 bytes
Desc: not available
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20150503/f133b272/attachment-0001.c>
-------------- next part --------------
A non-text attachment was scrubbed...
Name: makefile
Type: application/octet-stream
Size: 171 bytes
Desc: not available
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20150503/f133b272/attachment-0001.obj>


More information about the petsc-users mailing list