File: TSSetIJacobian.html

package info (click to toggle)
petsc 3.4.2.dfsg1-8.1
  • links: PTS, VCS
  • area: main
  • in suites: jessie, jessie-kfreebsd
  • size: 129,104 kB
  • ctags: 516,422
  • sloc: ansic: 395,939; cpp: 47,201; python: 34,788; makefile: 17,193; fortran: 16,251; f90: 1,592; objc: 954; sh: 822; xml: 621; java: 381; lisp: 293; csh: 241
file content (92 lines) | stat: -rw-r--r-- 6,104 bytes parent folder | download
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
<!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.4.2 2013-07-02</b></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 you 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,<A HREF="../Mat/MatStructure.html#MatStructure">MatStructure</A> *flag,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>flag </B></TD><TD>- flag indicating information about the preconditioner matrix
structure (same as flag in <A HREF="../KSP/KSPSetOperators.html#KSPSetOperators">KSPSetOperators</A>())
</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>
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>

<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>()
<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/ex22.c.html">src/ts/examples/tutorials/ex22.c.html</A><BR>
<A HREF="../../../src/ts/examples/tutorials/ex23.c.html">src/ts/examples/tutorials/ex23.c.html</A><BR>
<A HREF="../../../src/ts/examples/tutorials/ex24.c.html">src/ts/examples/tutorials/ex24.c.html</A><BR>
</BODY></HTML>