1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98
|
<!DOCTYPE HTML PUBLIC "-//W3C//DTD HTML3.2 EN">
<HTML>
<HEAD> <link rel="canonical" href="http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/TS/TSSetIJacobian.html" />
<META NAME="GENERATOR" CONTENT="DOCTEXT">
<TITLE>TSSetIJacobian</TITLE>
</HEAD>
<BODY BGCOLOR="FFFFFF">
<div id="version" align=right><b>petsc-3.7.5 2017-01-01</b></div>
<div id="bugreport" align=right><a href="mailto:petsc-maint@mcs.anl.gov?subject=Typo or Error in Documentation &body=Please describe the typo or error in the documentation: petsc-3.7.5 v3.7.5 docs/manualpages/TS/TSSetIJacobian.html "><small>Report Typos and Errors</small></a></div>
<A NAME="TSSetIJacobian"><H1>TSSetIJacobian</H1></A>
Set the function to compute the matrix dF/dU + a*dF/dU_t where F(t,U,U_t) is the function provided with <A HREF="../TS/TSSetIFunction.html#TSSetIFunction">TSSetIFunction</A>().
<H3><FONT COLOR="#CC3333">Synopsis</FONT></H3>
<PRE>
#include "petscts.h"
PetscErrorCode TSSetIJacobian(TS ts,Mat Amat,Mat Pmat,TSIJacobian f,void *ctx)
</PRE>
Logically Collective on <A HREF="../TS/TS.html#TS">TS</A>
<P>
<H3><FONT COLOR="#CC3333">Input Parameters</FONT></H3>
<TABLE border="0" cellpadding="0" cellspacing="0">
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>ts </B></TD><TD>- the <A HREF="../TS/TS.html#TS">TS</A> context obtained from <A HREF="../TS/TSCreate.html#TSCreate">TSCreate</A>()
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>Amat </B></TD><TD>- (approximate) Jacobian matrix
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>Pmat </B></TD><TD>- matrix used to compute preconditioner (usually the same as Amat)
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>f </B></TD><TD>- the Jacobian evaluation routine
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>ctx </B></TD><TD>- user-defined context for private data for the Jacobian evaluation routine (may be NULL)
</TD></TR></TABLE>
<P>
<H3><FONT COLOR="#CC3333">Calling sequence of f</FONT></H3>
<pre>
f(<A HREF="../TS/TS.html#TS">TS</A> ts,<A HREF="../Sys/PetscReal.html#PetscReal">PetscReal</A> t,<A HREF="../Vec/Vec.html#Vec">Vec</A> U,<A HREF="../Vec/Vec.html#Vec">Vec</A> U_t,<A HREF="../Sys/PetscReal.html#PetscReal">PetscReal</A> a,<A HREF="../Mat/Mat.html#Mat">Mat</A> Amat,<A HREF="../Mat/Mat.html#Mat">Mat</A> Pmat,void *ctx);
</pre>
<P>
<TABLE border="0" cellpadding="0" cellspacing="0">
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>t </B></TD><TD>- time at step/stage being solved
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>U </B></TD><TD>- state vector
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>U_t </B></TD><TD>- time derivative of state vector
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>a </B></TD><TD>- shift
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>Amat </B></TD><TD>- (approximate) Jacobian of F(t,U,W+a*U), equivalent to dF/dU + a*dF/dU_t
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>Pmat </B></TD><TD>- matrix used for constructing preconditioner, usually the same as Amat
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>ctx </B></TD><TD>- [optional] user-defined context for matrix evaluation routine
</TD></TR></TABLE>
<P>
<H3><FONT COLOR="#CC3333">Notes</FONT></H3>
The matrices Amat and Pmat are exactly the matrices that are used by <A HREF="../SNES/SNES.html#SNES">SNES</A> for the nonlinear solve.
<P>
If you know the operator Amat has a null space you can use <A HREF="../Mat/MatSetNullSpace.html#MatSetNullSpace">MatSetNullSpace</A>() and <A HREF="../Mat/MatSetTransposeNullSpace.html#MatSetTransposeNullSpace">MatSetTransposeNullSpace</A>() to supply the null
space to Amat and the <A HREF="../KSP/KSP.html#KSP">KSP</A> solvers will automatically use that null space as needed during the solution process.
<P>
The matrix dF/dU + a*dF/dU_t you provide turns out to be
the Jacobian of F(t,U,W+a*U) where F(t,U,U_t) = 0 is the DAE to be solved.
The time integrator internally approximates U_t by W+a*U where the positive "shift"
a and vector W depend on the integration method, step size, and past states. For example with
the backward Euler method a = 1/dt and W = -a*U(previous timestep) so
W + a*U = a*(U - U(previous timestep)) = (U - U(previous timestep))/dt
<P>
You must set all the diagonal entries of the matrices, if they are zero you must still set them with a zero value
<P>
The <A HREF="../TS/TS.html#TS">TS</A> solver may modify the nonzero structure and the entries of the matrices Amat and Pmat between the calls to f()
You should not assume the values are the same in the next call to f() as you set them in the previous call.
<P>
<P>
<H3><FONT COLOR="#CC3333">Keywords</FONT></H3>
<A HREF="../TS/TS.html#TS">TS</A>, timestep, DAE, Jacobian
<BR>
<P>
<H3><FONT COLOR="#CC3333">See Also</FONT></H3>
<A HREF="../TS/TSSetIFunction.html#TSSetIFunction">TSSetIFunction</A>(), <A HREF="../TS/TSSetRHSJacobian.html#TSSetRHSJacobian">TSSetRHSJacobian</A>(), <A HREF="../SNES/SNESComputeJacobianDefaultColor.html#SNESComputeJacobianDefaultColor">SNESComputeJacobianDefaultColor</A>(), <A HREF="../SNES/SNESComputeJacobianDefault.html#SNESComputeJacobianDefault">SNESComputeJacobianDefault</A>(), <A HREF="../TS/TSSetRHSFunction.html#TSSetRHSFunction">TSSetRHSFunction</A>()
<BR>
<P>
<P><B><P><B><FONT COLOR="#CC3333">Level:</FONT></B>beginner
<BR><FONT COLOR="#CC3333">Location:</FONT></B><A HREF="../../../src/ts/interface/ts.c.html#TSSetIJacobian">src/ts/interface/ts.c</A>
<BR><A HREF="./index.html">Index of all TS routines</A>
<BR><A HREF="../../index.html">Table of Contents for all manual pages</A>
<BR><A HREF="../singleindex.html">Index of all manual pages</A>
<P><H3><FONT COLOR="#CC3333">Examples</FONT></H3>
<A HREF="../../../src/ts/examples/tutorials/ex6.c.html">src/ts/examples/tutorials/ex6.c.html</A><BR>
<A HREF="../../../src/ts/examples/tutorials/ex8.c.html">src/ts/examples/tutorials/ex8.c.html</A><BR>
<A HREF="../../../src/ts/examples/tutorials/ex9.c.html">src/ts/examples/tutorials/ex9.c.html</A><BR>
<A HREF="../../../src/ts/examples/tutorials/ex10.c.html">src/ts/examples/tutorials/ex10.c.html</A><BR>
<A HREF="../../../src/ts/examples/tutorials/ex14.c.html">src/ts/examples/tutorials/ex14.c.html</A><BR>
<A HREF="../../../src/ts/examples/tutorials/ex15.c.html">src/ts/examples/tutorials/ex15.c.html</A><BR>
<A HREF="../../../src/ts/examples/tutorials/ex16.c.html">src/ts/examples/tutorials/ex16.c.html</A><BR>
<A HREF="../../../src/ts/examples/tutorials/ex17.c.html">src/ts/examples/tutorials/ex17.c.html</A><BR>
<A HREF="../../../src/ts/examples/tutorials/ex19.c.html">src/ts/examples/tutorials/ex19.c.html</A><BR>
<A HREF="../../../src/ts/examples/tutorials/ex20.c.html">src/ts/examples/tutorials/ex20.c.html</A><BR>
<A HREF="../../../src/ts/examples/tutorials/ex22.c.html">src/ts/examples/tutorials/ex22.c.html</A><BR>
</BODY></HTML>
|