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
|
<!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/TSSetCostIntegrand.html" />
<META NAME="GENERATOR" CONTENT="DOCTEXT">
<TITLE>TSSetCostIntegrand</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/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"
PetscErrorCode TSSetCostIntegrand(TS ts,PetscInt numcost,PetscErrorCode (*rf)(TS,PetscReal,Vec,Vec,void*),
PetscErrorCode (*drdyf)(TS,PetscReal,Vec,Vec*,void*),
PetscErrorCode (*drdpf)(TS,PetscReal,Vec,Vec*,void*),
PetscBool 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>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,NULL if not a function 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, NULL if not a function of p
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>fwd - flag indicating whether to evaluate cost integral in the forward run or the adjoint run</B></TD><TD>- - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
</TD></TR>
<P>
<H3><FONT COLOR="#CC3333">Calling sequence of rf</FONT></H3>
<pre>
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>
<TABLE border="0" cellpadding="0" cellspacing="0">
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>t </B></TD><TD>- current timestep
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>y </B></TD><TD>- input vector
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>f </B></TD><TD>- function result; one vector entry for each cost function
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>ctx </B></TD><TD>- [optional] user-defined function context
</TD></TR></TABLE>
<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>
Notes: For optimization there is generally 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="../TS/TSAdjointSetRHSJacobian.html#TSAdjointSetRHSJacobian">TSAdjointSetRHSJacobian</A>(),<A HREF="../TS/TSGetCostGradients.html#TSGetCostGradients">TSGetCostGradients</A>(), <A HREF="../TS/TSSetCostGradients.html#TSSetCostGradients">TSSetCostGradients</A>()
<BR><P><B><P><B><FONT COLOR="#CC3333">Level:</FONT></B>intermediate
<BR><FONT COLOR="#CC3333">Location:</FONT></B><A HREF="../../../src/ts/interface/ts.c.html#TSSetCostIntegrand">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>
</BODY></HTML>
|