[petsc-users] Jacobian finite difference approximation using coloring
Barry Smith
bsmith at mcs.anl.gov
Wed Nov 2 14:23:08 CDT 2011
It should be
ierr = SNESSetJacobian(sys.snes, sys.J, sys.P, SNESDefaultComputeJacobianColor, fdcoloring);
not
ierr = SNESSetJacobian(sys.snes, sys.J, sys.P, SNESDefaultComputeJacobianColor, &fdcoloring);
Barry
This is kind of our fault for having no type checking on contexts so the compiler cannot flag the problem.
On Nov 2, 2011, at 12:37 PM, Konstantinos Kontzialis wrote:
> Dear all,
>
> I am trying to compute the boundary layer over a flat plate using the discontinous galerkin method.
>
> I use the following sequence for computing the jacobian of a system using matrix coloring:
>
> ISColoring iscoloring;
> MatFDColoring fdcoloring;
>
> ierr = jacobian_diff_numerical(sys, &sys.P); /* Here I initialize the nonzerostructure of the matrix*/
> CHKERRQ(ierr);
>
> ierr = MatGetColoring(sys.P, MATCOLORING_ID, &iscoloring);
> CHKERRQ(ierr);
>
> ierr = MatFDColoringCreate(sys.P, iscoloring, &fdcoloring);
> CHKERRQ(ierr);
>
> ierr = MatFDColoringSetFunction(fdcoloring, base_residual_implicit,
> &sys);
> CHKERRQ(ierr);
>
> ierr = SNESSetJacobian(sys.snes, sys.J, sys.P, SNESDefaultComputeJacobianColor, &fdcoloring);
> CHKERRQ(ierr);
>
> I run my code as follows:
>
> mpiexec -n 8 ./hoac blasius -llf_flux -n_out 1 -end_time 10000.0 -implicit -implicit_type 3 -pc_type asm -snes_mf_operator -snes_max_fail 500 -snes_monitor -snes_stol 1.0e-50 -ksp_right_pc -sub_pc_type ilu -snes_converged_reason -ksp_gmres_restart 30 -snes_max_linear_solve_fail 10 -snes_max_it 1000 -sub_pc_factor_mat_ordering_type rcm -dt 5.e-4 -snes_rtol 1.0e-8 -gl -snes_converged_reason -ksp_converged_reason -ksp_monitor_true_residual -ksp_rtol 1.0e-12 -snes_atol 1.0e-6 -snes_ls_maxstep 5
>
>
> and I get this:
>
> Timestep 0: dt = 0.0005, T = 0, Res[rho] = 5.5383e-17, Res[rhou] = 0.0177116, Res[rhov] = 8.06867e-06, Res[E] = 9.04882e-06, CFL = 99.9994
>
> /*********************Stage 1 of SSPIRK (3,4)******************/
>
> 0 SNES Function norm 3.861205145119e-01
> [6]PETSC ERROR: --------------------- Error Message ------------------------------------
> [0]PETSC ERROR: --------------------- Error Message ------------------------------------
> [0]PETSC ERROR: Invalid argument!
> [0]PETSC ERROR: Wrong type of object: Parameter # 1!
> [0]PETSC ERROR: ------------------------------------------------------------------------
> [0]PETSC ERROR: Petsc Release Version 3.1.0, Patch 8, Thu Mar 17 13:37:48 CDT 2011
> [0]PETSC ERROR: See docs/changes/index.html for recent updates.
> [0]PETSC ERROR: See docs/faq.html for hints about trouble shooting.
> [0]PETSC ERROR: See docs/index.html for manual pages.
> [0]PETSC ERROR: ------------------------------------------------------------------------
> [0]PETSC ERROR: ./hoac on a linux-gnu named PlusSodaL by kontzialis Wed Nov 2 19:32:57 2011
> [0]PETSC ERROR: Libraries linked from /home/kontzialis/petsc-3.1-p8/linux-gnu-c-debug/lib
> [0]PETSC ERROR: Configure run at Tue Sep 27 13:09:04 2011
> [0]PETSC ERROR: Configure options --with-debugging=1 --with-shared=1 --with-shared-libraries --with-large-file-io=1 --with-precision=double --with-blacs=1 --download-blacs=yes --download-f-blas-lapack=yes --with-plapack=1 --download-plapack=yes --with-scalapack=1 --download-scalapack=yes --with-superlu=1 --download-superlu=yes --with-superlu_dist=1 --download-superlu_dist=yes --with-ml=1 --download-ml=yes --with-umfpack=1 --download-umfpack=yes --with-sundials=1 --download-sundials=1 --with-parmetis=1 --download-parmetis=1 --with-hypre=1 --download-hypre=1 --with-mpi-dir=/usr/lib/mpich2/bin
> [0]PETSC ERROR: ------------------------------------------------------------------------
> [0]PETSC ERROR: MatFDColoringGetFunction() line 205 in src/mat/matfd/fdmatrix.c
> [0]PETSC ERROR: SNESDefaultComputeJacobianColor() line 44 in src/snes/interface/snesj2.c
> [0]PETSC ERROR: SNESComputeJacobian() line 1198 in src/snes/interface/snes.c
> [0]PETSC ERROR: SNESSolve_LS() line 189 in src/snes/impls/ls/ls.c
> [0]PETSC ERROR: SNESSolve() line 2255 in src/snes/interface/snes.c
> [0]PETSC ERROR: User provided function() line 65 in "unknowndirectory/"../src/implicit_solve.c
> [0]PETSC ERROR: User provided function() line 215 in "unknowndirectory/"../src/implicit_time.c
> [0]PETSC ERROR: User provided function() line 1260 in "unknowndirectory/"../src/hoac.c
>
>
> What am I doing wrong here?
>
> Thank you,
>
> Kostas
More information about the petsc-users
mailing list