File: VecBoundGradientProjection.html

package info (click to toggle)
petsc 3.7.5%2Bdfsg1-4
  • links: PTS, VCS
  • area: main
  • in suites: stretch
  • size: 163,864 kB
  • ctags: 618,438
  • sloc: ansic: 515,133; python: 29,793; makefile: 20,458; fortran: 18,998; cpp: 6,515; f90: 3,914; sh: 1,012; xml: 621; objc: 445; csh: 240; java: 13
file content (109 lines) | stat: -rw-r--r-- 6,065 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
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
<!DOCTYPE HTML PUBLIC "-//W3C//DTD HTML3.2 EN">
<HTML>
<HEAD> <link rel="canonical" href="http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecBoundGradientProjection.html" />
<META NAME="GENERATOR" CONTENT="DOCTEXT">
<TITLE>VecBoundGradientProjection</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/Vec/VecBoundGradientProjection.html "><small>Report Typos and Errors</small></a></div>
<A NAME="VecBoundGradientProjection"><H1>VecBoundGradientProjection</H1></A>
Projects  vector according to this definition. If XL[i] &lt; X[i] &lt; XU[i], then GP[i] = G[i]; If X[i]&lt;=XL[i], then GP[i] = min(G[i],0); If X[i]&gt;=XU[i], then GP[i] = max(G[i],0); 
<H3><FONT COLOR="#CC3333">Synopsis</FONT></H3>
<PRE>
#include "petscvec.h"  
PetscErrorCode VecBoundGradientProjection(Vec G, Vec X, Vec XL, Vec XU, Vec GP)
</PRE>
<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>G </B></TD><TD>- current gradient vector
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>X </B></TD><TD>- current solution vector
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>XL </B></TD><TD>- lower bounds
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>XU </B></TD><TD>- upper bounds
</TD></TR></TABLE>
<P>
<H3><FONT COLOR="#CC3333">Output Parameter</FONT></H3>
<DT><B>GP </B> -gradient projection vector
<br>
<P>

C@*/
<A HREF="../Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</A> <A HREF="../Vec/VecBoundGradientProjection.html#VecBoundGradientProjection">VecBoundGradientProjection</A>(<A HREF="../Vec/Vec.html#Vec">Vec</A> G, <A HREF="../Vec/Vec.html#Vec">Vec</A> X, <A HREF="../Vec/Vec.html#Vec">Vec</A> XL, <A HREF="../Vec/Vec.html#Vec">Vec</A> XU, <A HREF="../Vec/Vec.html#Vec">Vec</A> GP)
{
<P>
<A HREF="../Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</A> ierr;
<A HREF="../Sys/PetscInt.html#PetscInt">PetscInt</A>       n,i;
<A HREF="../Sys/PetscReal.html#PetscReal">PetscReal</A>      *xptr,*xlptr,*xuptr,*gptr,*gpptr;
<A HREF="../Sys/PetscReal.html#PetscReal">PetscReal</A>      xval,gpval;
<P>
/* Project variables at the lower and upper bound */
<A HREF="../Sys/PetscFunctionBegin.html#PetscFunctionBegin">PetscFunctionBegin</A>;
PetscValidHeaderSpecific(G,VEC_CLASSID,1);
PetscValidHeaderSpecific(X,VEC_CLASSID,2);
PetscValidHeaderSpecific(XL,VEC_CLASSID,3);
PetscValidHeaderSpecific(XU,VEC_CLASSID,4);
PetscValidHeaderSpecific(GP,VEC_CLASSID,5);
<P>
ierr = <A HREF="../Vec/VecGetLocalSize.html#VecGetLocalSize">VecGetLocalSize</A>(X,&amp;n);<A HREF="../Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</A>(ierr);
<P>
ierr=<A HREF="../Vec/VecGetArray.html#VecGetArray">VecGetArray</A>(X,&amp;xptr);<A HREF="../Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</A>(ierr);
ierr=<A HREF="../Vec/VecGetArray.html#VecGetArray">VecGetArray</A>(XL,&amp;xlptr);<A HREF="../Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</A>(ierr);
ierr=<A HREF="../Vec/VecGetArray.html#VecGetArray">VecGetArray</A>(XU,&amp;xuptr);<A HREF="../Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</A>(ierr);
ierr=<A HREF="../Vec/VecGetArray.html#VecGetArray">VecGetArray</A>(G,&amp;gptr);<A HREF="../Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</A>(ierr);
if (G!=GP){
ierr=<A HREF="../Vec/VecGetArray.html#VecGetArray">VecGetArray</A>(GP,&amp;gpptr);<A HREF="../Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</A>(ierr);
} else { gpptr=gptr; }
<P>
for (i=0; i&lt;n; ++i){
gpval = gptr[i]; xval = xptr[i];
<P>
if (gpval&gt;0 &amp;&amp; xval&lt;=xlptr[i]){
gpval = 0;
} else if (gpval&lt;0 &amp;&amp; xval&gt;=xuptr[i]){
gpval = 0;
}
gpptr[i] = gpval;
}
<P>
ierr=<A HREF="../Vec/VecRestoreArray.html#VecRestoreArray">VecRestoreArray</A>(X,&amp;xptr);<A HREF="../Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</A>(ierr);
ierr=<A HREF="../Vec/VecRestoreArray.html#VecRestoreArray">VecRestoreArray</A>(XL,&amp;xlptr);<A HREF="../Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</A>(ierr);
ierr=<A HREF="../Vec/VecRestoreArray.html#VecRestoreArray">VecRestoreArray</A>(XU,&amp;xuptr);<A HREF="../Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</A>(ierr);
ierr=<A HREF="../Vec/VecRestoreArray.html#VecRestoreArray">VecRestoreArray</A>(G,&amp;gptr);<A HREF="../Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</A>(ierr);
if (G!=GP){
ierr=<A HREF="../Vec/VecRestoreArray.html#VecRestoreArray">VecRestoreArray</A>(GP,&amp;gpptr);<A HREF="../Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</A>(ierr);
}
<A HREF="../Sys/PetscFunctionReturn.html#PetscFunctionReturn">PetscFunctionReturn</A>(0);
}
#endif
<P>
#undef __FUNCT__
#define __FUNCT__ "VecStepMaxBounded"
/*@
VecStepMaxBounded - See below
<P>
Collective on <A HREF="../Vec/Vec.html#Vec">Vec</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>X  </B></TD><TD>- vector with no negative entries
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>XL </B></TD><TD>- lower bounds
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>XU </B></TD><TD>- upper bounds
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>DX  </B></TD><TD>- step direction, can have negative, positive or zero entries
</TD></TR></TABLE>
<P>
<H3><FONT COLOR="#CC3333">Output Parameter</FONT></H3>
<DT><B>stepmax </B> -minimum value so that X[i] + stepmax*DX[i] &lt;= XL[i]  or  XU[i] &lt;= X[i] + stepmax*DX[i]
<br>
<P>
<P><B><P><B><FONT COLOR="#CC3333">Level:</FONT></B>advanced
<BR><FONT COLOR="#CC3333">Location:</FONT></B><A HREF="../../../src/vec/vec/utils/projection.c.html#VecBoundGradientProjection">src/vec/vec/utils/projection.c</A>
<BR><A HREF="./index.html">Index of all Vec 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>