[petsc-users] Question about using MatSNESMFWPSetComputeNormU

Barry Smith bsmith at mcs.anl.gov
Thu Aug 18 14:27:38 CDT 2016


   You can't use that Jacobian function SNESComputeJacobianDefault  with matrix free, it tries to compute the matrix entries and stick them into the matrix. You can use MatMFFDComputeJacobian

> On Aug 18, 2016, at 2:03 PM, 최경준 <kyungjun.choi92 at gmail.com> wrote:
> 
> I got stuck at FormJacobian stage.
> 
> - call SNESComputeJacobianDefault(snes, v, J, pJ, FormResidual, ier)  --> J & pJ are same with A matrix-free matrix (input argument)
> 
> 
> 
> with these kind of messages..
> 
> [0]PETSC ERROR: No support for this operation for this object type
> [0]PETSC ERROR: Mat type mffd
> 
> 
> 
> Guess it's because I used A matrix-free matrix (which is mffd type) into pJ position.
> 
> Is there any solution for this kind of situation?
> 
> 
> 2016-08-19 2:05 GMT+09:00 Matthew Knepley <knepley at gmail.com>:
> On Thu, Aug 18, 2016 at 12:04 PM, 최경준 <kyungjun.choi92 at gmail.com> wrote:
> Then in order not to use preconditioner, 
> 
> is it ok if I just put A matrix-free matrix (made from MatCreateSNESMF()) into the place where preA should be?
> 
> Yes, but again the solve will likely perform very poorly.
> 
>   Thanks,
> 
>      Matt
>  
> The flow goes like this
> - call SNESCreate
> - call SNESSetFunction(snes, r, FormResidual, userctx, ier)
> - call MatCreateSNESMF(snes, A, ier)
> - call SNESSetJacobian(snes, A, A, FormJacobian, userctx, ier)
> - call SNESSetFromOptions()
> 
> - call SNESGetKSP(snes, ksp, ier)
> - call KSPSetType(ksp, KSPGMRES, ier)
> - call KSPGetPC(ksp, pc, ier)
> - call PCSetType(pc, PCNONE, ier)
> - call KSPGMRESSetRestart(ksp, 30, ier)
> 
> - call SNESSolve()
> .
> .
> 
> 
> and inside the FormJacobian routine
> - call SNESComputeJacobian(snes, v, J, pJ, userctx, ier)  --> J and pJ must be pointed with A and A.
> 
> 
> 
> Thank you again,
> 
> Kyungjun.
> 
> 2016-08-19 1:44 GMT+09:00 Matthew Knepley <knepley at gmail.com>:
> On Thu, Aug 18, 2016 at 11:42 AM, 최경준 <kyungjun.choi92 at gmail.com> wrote:
> Thanks for your helpful answers.
> 
> Here's another question...
> 
> As I read some example PETSc codes, I noticed that there should be a preconditioning matrix (e.g. approx. jacobian matrix) when using MatCreateSNESMF().
> 
> I mean,
> after calling MatCreateSNESMF(snes, A, ier),
> there should be another matrix preA(preconditioning matrix) to use SNESSetJacobian(snes, A, preA, FormJacobian, ctx, ier).
> 
> 
> 1) Is there any way that I can use matrix-free method without making preconditioning matrix?
> 
> Don't use a preconditioner. As you might expect, this does not often work out well.
>  
> 2) I have a reference code, and the code adopts 
> 
> MatFDColoringCreate()
> and finally uses
> SNESComputeJacobianDefaultColor() at FormJacobian stage.
> 
> But I can't see the inside of the fdcolor and I'm curious of this mechanism. Can you explain this very briefly or tell me an example code that I can refer to. ( I think none of PETSc example code is using fdcolor..)
> 
> This is the default, so there is no need for all that code. We use naive graph 2-coloring. I think there might be a review article by Alex Pothen about that.
> 
>   Thanks,
> 
>     Matt
>  
> 
> Best,
> 
> Kyungjun.
> 
> 2016-08-19 0:54 GMT+09:00 Matthew Knepley <knepley at gmail.com>:
> On Thu, Aug 18, 2016 at 10:39 AM, 최경준 <kyungjun.choi92 at gmail.com> wrote:
> 1) I wanna know the difference between applying option with command line and within source code.
> From my experience, command line option helps set other default settings that I didn't applied, I guess.
> 
> The command line arguments are applied to an object when *SetFromOptions() is called, so in this case
> you want SNESSetFromOptions() on the solver. There should be no difference from using the API.
>  
> 2) I made a matrix-free matrix with MatCreateSNESMF function, and every time I check my snes context with SNESView, 
> 
>     Mat Object:     1 MPI processes
>       type: mffd
>       rows=11616, cols=11616
>         Matrix-free approximation:
>           err=1.49012e-08 (relative error in function evaluation)
>           The compute h routine has not yet been set
> 
> at the end of line shows there's no routine for computing h value.
> I used MatMFFDWPSetComputeNormU function, but it didn't work I think.
> Is it ok if I leave the h value that way? Or should I have to set h computing routine?
> 
> I am guessing you are calling the function on a different object from the one that is viewed here.
> However, there will always be a default function for computing h.
> 
>   Thanks,
> 
>     Matt
>  
> Kyungjun.
> 
> 2016-08-18 23:18 GMT+09:00 Matthew Knepley <knepley at gmail.com>:
> On Thu, Aug 18, 2016 at 8:35 AM, 최경준 <kyungjun.choi92 at gmail.com> wrote:
> Hi, I'm trying to set my SNES matrix-free with Walker & Pernice way of computing h value.
> 
> I found above command (MatSNESMFWPSetComputeNormU) but my fortran compiler couldn't fine any reference of that command.
> 
> I checked Petsc changes log, but there weren't any mentions about that command.
> 
> Should I have to include another specific header file? 
> 
> We have this function
> 
>   http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatMFFDWPSetComputeNormU.html
> 
> but I would recommend using the command line option
> 
>   -mat_mffd_compute_normu
>  
>   Thanks,
> 
>      Matt
> 
> Thank you always.
> 
> 
> 
> -- 
> 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
> 
> 
> 
> 
> -- 
> 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
> 
> 
> 
> 
> -- 
> 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
> 
> 
> 
> 
> -- 
> 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
> 



More information about the petsc-users mailing list