File: PetscDSAddBoundary.html

package info (click to toggle)
petsc 3.14.5%2Bdfsg1-4
  • links: PTS, VCS
  • area: main
  • in suites: bullseye
  • size: 266,472 kB
  • sloc: ansic: 680,898; python: 33,303; cpp: 16,324; makefile: 14,022; f90: 13,731; javascript: 10,713; fortran: 9,581; sh: 1,373; xml: 619; objc: 445; csh: 192; pascal: 148; java: 13
file content (123 lines) | stat: -rw-r--r-- 9,363 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
110
111
112
113
114
115
116
117
118
119
120
121
122
123
<!DOCTYPE HTML PUBLIC "-//W3C//DTD HTML3.2 EN">
<HTML>
<HEAD> <link rel="canonical" href="http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DT/PetscDSAddBoundary.html" />
<META NAME="GENERATOR" CONTENT="DOCTEXT">
<TITLE>PetscDSAddBoundary</TITLE>
</HEAD>
<BODY BGCOLOR="FFFFFF">
   <div id="version" align=right><b>petsc-3.14.5 2021-03-03</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.14.5 v3.14.5 docs/manualpages/DT/PetscDSAddBoundary.html "><small>Report Typos and Errors</small></a></div>
<A NAME="PetscDSAddBoundary"><H1>PetscDSAddBoundary</H1></A>
Add a boundary condition to the model. The pointwise functions are used to provide boundary values for essential boundary conditions. In FEM, they are acting upon by dual basis functionals to generate FEM coefficients which are fixed. Natural boundary conditions signal to PETSc that boundary integrals should be performaed, using the kernels from <A HREF="../DT/PetscDSSetBdResidual.html#PetscDSSetBdResidual">PetscDSSetBdResidual</A>(). 
<H3><FONT COLOR="#CC3333">Synopsis</FONT></H3>
<PRE>
#include "petscds.h" 
<A HREF="../Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</A> <A HREF="../DT/PetscDSAddBoundary.html#PetscDSAddBoundary">PetscDSAddBoundary</A>(<A HREF="../DM/PetscDS.html#PetscDS">PetscDS</A> ds, <A HREF="../DM/DMBoundaryConditionType.html#DMBoundaryConditionType">DMBoundaryConditionType</A> type, const char name[], const char labelname[], <A HREF="../Sys/PetscInt.html#PetscInt">PetscInt</A> field, <A HREF="../Sys/PetscInt.html#PetscInt">PetscInt</A> numcomps, const <A HREF="../Sys/PetscInt.html#PetscInt">PetscInt</A> *comps, void (*bcFunc)(void), void (*bcFunc_t)(void), <A HREF="../Sys/PetscInt.html#PetscInt">PetscInt</A> numids, const <A HREF="../Sys/PetscInt.html#PetscInt">PetscInt</A> *ids, void *ctx)
</PRE>
Collective on ds
<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>ds          </B></TD><TD>- The <A HREF="../DM/PetscDS.html#PetscDS">PetscDS</A> object
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>type        </B></TD><TD>- The type of condition, e.g. <A HREF="../DM/DMBoundaryConditionType.html#DMBoundaryConditionType">DM_BC_ESSENTIAL</A>/<A HREF="../DM/DMBoundaryConditionType.html#DMBoundaryConditionType">DM_BC_ESSENTIAL_FIELD</A> (Dirichlet), or <A HREF="../DM/DMBoundaryConditionType.html#DMBoundaryConditionType">DM_BC_NATURAL</A> (Neumann)
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>name        </B></TD><TD>- The BC name
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>labelname   </B></TD><TD>- The label defining constrained points
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>field       </B></TD><TD>- The field to constrain
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>numcomps    </B></TD><TD>- The number of constrained field components (0 will constrain all fields)
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>comps       </B></TD><TD>- An array of constrained component numbers
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>bcFunc      </B></TD><TD>- A pointwise function giving boundary values
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>bcFunc_t    </B></TD><TD>- A pointwise function giving the time derviative of the boundary values, or NULL
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>numids      </B></TD><TD>- The number of <A HREF="../DM/DMLabel.html#DMLabel">DMLabel</A> ids for constrained points
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>ids         </B></TD><TD>- An array of ids for constrained points
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>ctx         </B></TD><TD>- An optional user context for bcFunc
</TD></TR></TABLE>
<P>
<H3><FONT COLOR="#CC3333">Options Database Keys</FONT></H3>
<TABLE border="0" cellpadding="0" cellspacing="0">
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-bc_&lt;boundary name&gt; &lt;num&gt; </B></TD><TD>- Overrides the boundary ids
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-bc_&lt;boundary name&gt;_comp &lt;num&gt; </B></TD><TD>- Overrides the boundary components
</TD></TR></TABLE>
<P>
<H3><FONT COLOR="#CC3333">Note</FONT></H3>
<H3><FONT COLOR="#CC3333">Both bcFunc abd bcFunc_t will depend on the boundary condition type. If the type if DM_BC_ESSENTIAL, Then the calling sequence is</FONT></H3>
<P>
<pre>
bcFunc(<A HREF="../Sys/PetscInt.html#PetscInt">PetscInt</A> dim, <A HREF="../Sys/PetscReal.html#PetscReal">PetscReal</A> time, const <A HREF="../Sys/PetscReal.html#PetscReal">PetscReal</A> x[], <A HREF="../Sys/PetscInt.html#PetscInt">PetscInt</A> Nc, <A HREF="../Sys/PetscScalar.html#PetscScalar">PetscScalar</A> bcval[])
</pre>
<P>
<H3><FONT COLOR="#CC3333">If the type is DM_BC_ESSENTIAL_FIELD or other _FIELD value, then the calling sequence is</FONT></H3>
<P>
<pre>
bcFunc(<A HREF="../Sys/PetscInt.html#PetscInt">PetscInt</A> dim, <A HREF="../Sys/PetscInt.html#PetscInt">PetscInt</A> Nf, <A HREF="../Sys/PetscInt.html#PetscInt">PetscInt</A> NfAux,
</pre>
<pre>
       const <A HREF="../Sys/PetscInt.html#PetscInt">PetscInt</A> uOff[], const <A HREF="../Sys/PetscInt.html#PetscInt">PetscInt</A> uOff_x[], const <A HREF="../Sys/PetscScalar.html#PetscScalar">PetscScalar</A> u[], const <A HREF="../Sys/PetscScalar.html#PetscScalar">PetscScalar</A> u_t[], const <A HREF="../Sys/PetscScalar.html#PetscScalar">PetscScalar</A> u_x[],
</pre>
<pre>
       const <A HREF="../Sys/PetscInt.html#PetscInt">PetscInt</A> aOff[], const <A HREF="../Sys/PetscInt.html#PetscInt">PetscInt</A> aOff_x[], const <A HREF="../Sys/PetscScalar.html#PetscScalar">PetscScalar</A> a[], const <A HREF="../Sys/PetscScalar.html#PetscScalar">PetscScalar</A> a_t[], const <A HREF="../Sys/PetscScalar.html#PetscScalar">PetscScalar</A> a_x[],
</pre>
<pre>
       <A HREF="../Sys/PetscReal.html#PetscReal">PetscReal</A> time, const <A HREF="../Sys/PetscReal.html#PetscReal">PetscReal</A> x[], <A HREF="../Sys/PetscScalar.html#PetscScalar">PetscScalar</A> bcval[])
</pre>
<P>
<TABLE border="0" cellpadding="0" cellspacing="0">
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>dim </B></TD><TD>- the spatial dimension
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>Nf </B></TD><TD>- the number of fields
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>uOff </B></TD><TD>- the offset into u[] and u_t[] for each field
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>uOff_x </B></TD><TD>- the offset into u_x[] for each field
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>u </B></TD><TD>- each field evaluated at the current point
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>u_t </B></TD><TD>- the time derivative of each field evaluated at the current point
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>u_x </B></TD><TD>- the gradient of each field evaluated at the current point
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>aOff </B></TD><TD>- the offset into a[] and a_t[] for each auxiliary field
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>aOff_x </B></TD><TD>- the offset into a_x[] for each auxiliary field
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>a </B></TD><TD>- each auxiliary field evaluated at the current point
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>a_t </B></TD><TD>- the time derivative of each auxiliary field evaluated at the current point
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>a_x </B></TD><TD>- the gradient of auxiliary each field evaluated at the current point
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>t </B></TD><TD>- current time
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>x </B></TD><TD>- coordinates of the current point
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>numConstants </B></TD><TD>- number of constant parameters
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>constants </B></TD><TD>- constant parameters
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>bcval </B></TD><TD>- output values at the current point
</TD></TR></TABLE>
<P>

<P>
<H3><FONT COLOR="#CC3333">See Also</FONT></H3>
 <A HREF="../DT/PetscDSGetBoundary.html#PetscDSGetBoundary">PetscDSGetBoundary</A>(), <A HREF="../DT/PetscDSSetResidual.html#PetscDSSetResidual">PetscDSSetResidual</A>(), <A HREF="../DT/PetscDSSetBdResidual.html#PetscDSSetBdResidual">PetscDSSetBdResidual</A>()
<BR><P><B></B><H3><FONT COLOR="#CC3333">Level</FONT></H3>developer<BR>
<H3><FONT COLOR="#CC3333">Location</FONT></H3>
</B><A HREF="../../../src/dm/dt/interface/dtds.c.html#PetscDSAddBoundary">src/dm/dt/interface/dtds.c</A>
<P><H3><FONT COLOR="#CC3333">Examples</FONT></H3>
<A HREF="../../../src/ts/tutorials/ex11.c.html">src/ts/tutorials/ex11.c.html</A><BR>
<BR><A HREF="./index.html">Index of all DT 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>