Please send a run with optimization turned on (--with-debugging=0 in 
./configure) and -log_view  without the actual timing information we are just 
guessing where the time is spent.

  If your problem has a natural block size then using baij should be a bit 
faster than aij, but not dramatically better

   Barry


> On Dec 2, 2019, at 4:30 AM, Li Luo <li....@kaust.edu.sa> wrote:
> 
> Thank you very much! It looks forming an analytic Jacobian is the only choice.
> 
> Best,
> Li
> 
> On Mon, Dec 2, 2019 at 3:21 PM Matthew Knepley <knep...@gmail.com> wrote:
> On Mon, Dec 2, 2019 at 4:04 AM Li Luo <li....@kaust.edu.sa> wrote:
> Thank you for your reply.
> 
> The matrix is small with only 67500 rows, but is relatively dense since a 
> second-order discontinuous Galerkin FEM is used, nonzeros=23,036,400.
> 
> This is very dense, 0.5% fill or 340 nonzeros per row.
>  
> The number of colors is 539 as shown by using -mat_fd_coloring_view:
> 
> Coloring is not appropriate for this matrix since you have enormous dense 
> blocks (I am guessing). It could work if you statically
> condense them out or had a fast analytic Jacobian. With 540 colors, it takes 
> 540 matvecs to generate the action of the Jacobian.
> 
>   Thanks,
> 
>      Matt
>  
> MatFDColoring Object: 64 MPI processes
>   type not yet set
>   Error tolerance=1.49012e-08
>   Umin=1.49012e-06
>   Number of colors=539
>   Information for color 0
>     Number of columns 1
>       378
>     Number of rows 756
>       0 1188
>       1 1188
>       2 1188
>       3 1188
>       4 1188
>       5 1188
>       ...
> 
> Is this normal?
> When using MCFD, is there any difference using mpiaij and mpibaij?
> 
> Best,
> Li
> 
> On Mon, Dec 2, 2019 at 10:03 AM Smith, Barry F. <bsm...@mcs.anl.gov> wrote:
> 
>   How many colors is it requiring?   And how long is the MatGetColoring() 
> taking? Are you running in parallel?  The MatGetColoring() MATCOLORINGSL uses 
> a sequential coloring algorithm so if your matrix is large and parallel the 
> coloring will take a long time. The parallel colorings are MATCOLORINGGREEDY 
> and MATCOLORINGJP
> 
>   Barry
> 
> 
> > On Dec 1, 2019, at 12:56 AM, Li Luo <li....@kaust.edu.sa> wrote:
> > 
> > Dear Developers,
> > 
> > I tried to use the multi-color finite-difference (MC-FD) method for 
> > constructing the Jacobians. However, I find it is very slow compared to the 
> > exact Jacobian. 
> > My implementation of MC-FD Jacobian is posted below, would you please check 
> > whether I am correct? Anything missed? Thank you!
> > 
> > ////////// Setup phase:
> >           MatStructure flag;
> >           ISColoring   iscoloring;
> >           ierr = MatGetColoring(Jac,MATCOLORINGSL,&iscoloring);
> >           ierr = MatFDColoringCreate(Jac,iscoloring,&this->matfdcoloring);
> >           ierr = 
> > MatFDColoringSetFunction(this->matfdcoloring,(PetscErrorCode 
> > (*)(void))__libmesh_petsc_snes_residual,(void *)this);
> >           ierr = MatFDColoringSetFromOptions(this->matfdcoloring);
> >           ierr = ISColoringDestroy(&iscoloring);
> > 
> > //////////// Apply:
> >           ierr = MatZeroEntries(*jac);CHKERRQ(ierr);
> >           ierr = 
> > MatFDColoringApply(*jac,solver->matfdcoloring,x,msflag,snes);
> > 
> > Best regards,
> > Li Luo
> > 
> > This message and its contents, including attachments are intended solely 
> > for the original recipient. If you are not the intended recipient or have 
> > received this message in error, please notify me immediately and delete 
> > this message from your computer system. Any unauthorized use or 
> > distribution is prohibited. Please consider the environment before printing 
> > this email.
> 
> 
> 
> -- 
> Postdoctoral Fellow
> Extreme Computing Research Center
> King Abdullah University of Science & Technology
> https://sites.google.com/site/rolyliluo/
> 
> This message and its contents, including attachments are intended solely for 
> the original recipient. If you are not the intended recipient or have 
> received this message in error, please notify me immediately and delete this 
> message from your computer system. Any unauthorized use or distribution is 
> prohibited. Please consider the environment before printing this email.
> 
> 
> -- 
> 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
> 
> https://www.cse.buffalo.edu/~knepley/
> 
> 
> -- 
> Postdoctoral Fellow
> Extreme Computing Research Center
> King Abdullah University of Science & Technology
> https://sites.google.com/site/rolyliluo/
> 
> This message and its contents, including attachments are intended solely for 
> the original recipient. If you are not the intended recipient or have 
> received this message in error, please notify me immediately and delete this 
> message from your computer system. Any unauthorized use or distribution is 
> prohibited. Please consider the environment before printing this email.

Reply via email to