File: TSSetCostIntegrand.html

package info (click to toggle)
petsc 3.10.3%2Bdfsg1-5
  • links: PTS, VCS
  • area: main
  • in suites: buster
  • size: 209,064 kB
  • sloc: ansic: 587,333; python: 29,696; makefile: 12,445; fortran: 11,626; f90: 9,677; cpp: 8,768; sh: 1,027; xml: 621; objc: 445; csh: 194; java: 13
file content (73 lines) | stat: -rw-r--r-- 5,758 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
<!DOCTYPE HTML PUBLIC "-//W3C//DTD HTML3.2 EN">
<HTML>
<HEAD> <link rel="canonical" href="http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Sensitivity/TSSetCostIntegrand.html" />
<META NAME="GENERATOR" CONTENT="DOCTEXT">
<TITLE>TSSetCostIntegrand</TITLE>
</HEAD>
<BODY BGCOLOR="FFFFFF">
   <div id="version" align=right><b>petsc-3.10.3 2018-12-18</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.10.3 v3.10.3 docs/manualpages/Sensitivity/TSSetCostIntegrand.html "><small>Report Typos and Errors</small></a></div>
<A NAME="TSSetCostIntegrand"><H1>TSSetCostIntegrand</H1></A>
Sets the routine for evaluating the integral term in one or more cost functions 
<H3><FONT COLOR="#CC3333">Synopsis</FONT></H3>
<PRE>
#include "petscts.h"  
<A HREF="../Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</A> <A HREF="../Sensitivity/TSSetCostIntegrand.html#TSSetCostIntegrand">TSSetCostIntegrand</A>(<A HREF="../TS/TS.html#TS">TS</A> ts,<A HREF="../Sys/PetscInt.html#PetscInt">PetscInt</A> numcost,<A HREF="../Vec/Vec.html#Vec">Vec</A> costintegral,<A HREF="../Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</A> (*rf)(<A HREF="../TS/TS.html#TS">TS</A>,<A HREF="../Sys/PetscReal.html#PetscReal">PetscReal</A>,<A HREF="../Vec/Vec.html#Vec">Vec</A>,<A HREF="../Vec/Vec.html#Vec">Vec</A>,void*),
                                                          <A HREF="../Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</A> (*drdyf)(<A HREF="../TS/TS.html#TS">TS</A>,<A HREF="../Sys/PetscReal.html#PetscReal">PetscReal</A>,<A HREF="../Vec/Vec.html#Vec">Vec</A>,<A HREF="../Vec/Vec.html#Vec">Vec</A>*,void*),
                                                          <A HREF="../Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</A> (*drdpf)(<A HREF="../TS/TS.html#TS">TS</A>,<A HREF="../Sys/PetscReal.html#PetscReal">PetscReal</A>,<A HREF="../Vec/Vec.html#Vec">Vec</A>,<A HREF="../Vec/Vec.html#Vec">Vec</A>*,void*),
                                                          <A HREF="../Sys/PetscBool.html#PetscBool">PetscBool</A> fwd,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>numcost </B></TD><TD>- number of gradients to be computed, this is the number of cost functions
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>costintegral </B></TD><TD>- vector that stores the integral values
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>rf </B></TD><TD>- routine for evaluating the integrand function
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>drdyf </B></TD><TD>- function that computes the gradients of the r's with respect to y
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>drdpf </B></TD><TD>- function that computes the gradients of the r's with respect to p, can be NULL if parametric sensitivity is not desired (mu=NULL)
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>fwd </B></TD><TD>- flag indicating whether to evaluate cost integral in the forward run or the adjoint run
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>ctx </B></TD><TD>- [optional] user-defined context for private data for the function evaluation routine (may be NULL)
</TD></TR></TABLE>
<P>
<H3><FONT COLOR="#CC3333">Calling sequence of rf</FONT></H3>
<pre>
  <A HREF="../Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</A> rf(<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> y,<A HREF="../Vec/Vec.html#Vec">Vec</A> f,void *ctx);
</pre>
<P>
<H3><FONT COLOR="#CC3333">Calling sequence of drdyf</FONT></H3>
<pre>
  PetscErroCode drdyf(<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> y,<A HREF="../Vec/Vec.html#Vec">Vec</A> *drdy,void *ctx);
</pre>
<P>
<H3><FONT COLOR="#CC3333">Calling sequence of drdpf</FONT></H3>
<pre>
  PetscErroCode drdpf(<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> y,<A HREF="../Vec/Vec.html#Vec">Vec</A> *drdp,void *ctx);
</pre>
<P>

<P>
<H3><FONT COLOR="#CC3333">Notes</FONT></H3>
For optimization there is usually a single cost function (numcost = 1). For sensitivities there may be multiple cost functions
<P>
<H3><FONT COLOR="#CC3333">Keywords</FONT></H3>
 <A HREF="../TS/TS.html#TS">TS</A>, sensitivity analysis, timestep, set, quadrature, function
<BR>
<P>
<H3><FONT COLOR="#CC3333">See Also</FONT></H3>
 <A HREF="../Sensitivity/TSSetRHSJacobianP.html#TSSetRHSJacobianP">TSSetRHSJacobianP</A>(), <A HREF="../Sensitivity/TSGetCostGradients.html#TSGetCostGradients">TSGetCostGradients</A>(), <A HREF="../Sensitivity/TSSetCostGradients.html#TSSetCostGradients">TSSetCostGradients</A>()
<BR><P><B></B><H3><FONT COLOR="#CC3333">Level</FONT></H3>intermediate<BR>
<H3><FONT COLOR="#CC3333">Location</FONT></H3>
</B><A HREF="../../../src/ts/interface/sensitivity/tssen.c.html#TSSetCostIntegrand">src/ts/interface/sensitivity/tssen.c</A>
<BR><A HREF="./index.html">Index of all Sensitivity 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>
</BODY></HTML>