File: PetscFEIntegrateJacobian.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 (178 lines) | stat: -rw-r--r-- 10,535 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
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
<!DOCTYPE HTML PUBLIC "-//W3C//DTD HTML3.2 EN">
<HTML>
<HEAD> <link rel="canonical" href="http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DM/PetscFEIntegrateJacobian.html" />
<META NAME="GENERATOR" CONTENT="DOCTEXT">
<TITLE>PetscFEIntegrateJacobian</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/DM/PetscFEIntegrateJacobian.html "><small>Report Typos and Errors</small></a></div>
<A NAME="PetscFEIntegrateJacobian"><H1>PetscFEIntegrateJacobian</H1></A>
Produce the element Jacobian for a chunk of elements by quadrature integration 
<H3><FONT COLOR="#CC3333">Synopsis</FONT></H3>
<PRE>
#include "petscfe.h" 
PetscErrorCode PetscFERefine(PetscFE fe, PetscFE *feRef)
</PRE>
Not collective
<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>fem          </B></TD><TD>- The <A HREF="../DM/PetscFE.html#PetscFE">PetscFE</A> object for the field being integrated
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>prob         </B></TD><TD>- The <A HREF="../DM/PetscDS.html#PetscDS">PetscDS</A> specifying the discretizations and continuum functions
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>jtype        </B></TD><TD>- The type of matrix pointwise functions that should be used
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>fieldI       </B></TD><TD>- The test field being integrated
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>fieldJ       </B></TD><TD>- The basis field being integrated
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>Ne           </B></TD><TD>- The number of elements in the chunk
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>geom         </B></TD><TD>- The cell geometry for each cell in the chunk
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>coefficients </B></TD><TD>- The array of FEM basis coefficients for the elements for the Jacobian evaluation point
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>coefficients_t </B></TD><TD>- The array of FEM basis time derivative coefficients for the elements
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>probAux      </B></TD><TD>- The <A HREF="../DM/PetscDS.html#PetscDS">PetscDS</A> specifying the auxiliary discretizations
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>coefficientsAux </B></TD><TD>- The array of FEM auxiliary basis coefficients for the elements
</TD></TR></TABLE>
<P>
Output Parameter
<DT><B>elemMat      </B> -the element matrices for the Jacobian from each element
<br>
<P>
<H3><FONT COLOR="#CC3333">Note</FONT></H3>
<pre>
Loop over batch of elements (e):
</pre>
<pre>
  Loop over element matrix entries (f,fc,g,gc --&gt; i,j):
</pre>
<pre>
    Loop over quadrature points (q):
</pre>
<pre>
      Make u_q and gradU_q (loops over fields,Nb,Ncomp)
</pre>
<pre>
        elemMat[i,j] += \psi^{fc}_f(q) g0_{fc,gc}(u, \nabla u) \phi^{gc}_g(q)
</pre>
<pre>
                     + \psi^{fc}_f(q) \cdot g1_{fc,gc,dg}(u, \nabla u) \nabla\phi^{gc}_g(q)
</pre>
<pre>
                     + \nabla\psi^{fc}_f(q) \cdot g2_{fc,gc,df}(u, \nabla u) \phi^{gc}_g(q)
</pre>
<pre>
                     + \nabla\psi^{fc}_f(q) \cdot g3_{fc,gc,df,dg}(u, \nabla u) \nabla\phi^{gc}_g(q)
</pre>
*/
<A HREF="../Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</A> <A HREF="../DM/PetscFEIntegrateJacobian.html#PetscFEIntegrateJacobian">PetscFEIntegrateJacobian</A>(<A HREF="../DM/PetscFE.html#PetscFE">PetscFE</A> fem, <A HREF="../DM/PetscDS.html#PetscDS">PetscDS</A> prob, <A HREF="../DM/PetscFEJacobianType.html#PetscFEJacobianType">PetscFEJacobianType</A> jtype, <A HREF="../Sys/PetscInt.html#PetscInt">PetscInt</A> fieldI, <A HREF="../Sys/PetscInt.html#PetscInt">PetscInt</A> fieldJ, <A HREF="../Sys/PetscInt.html#PetscInt">PetscInt</A> Ne, PetscFECellGeom *geom,
const <A HREF="../Sys/PetscScalar.html#PetscScalar">PetscScalar</A> coefficients[], const <A HREF="../Sys/PetscScalar.html#PetscScalar">PetscScalar</A> coefficients_t[], <A HREF="../DM/PetscDS.html#PetscDS">PetscDS</A> probAux, const <A HREF="../Sys/PetscScalar.html#PetscScalar">PetscScalar</A> coefficientsAux[], <A HREF="../Sys/PetscScalar.html#PetscScalar">PetscScalar</A> elemMat[])
{
<A HREF="../Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</A> ierr;
<P>
<A HREF="../Sys/PetscFunctionBegin.html#PetscFunctionBegin">PetscFunctionBegin</A>;
PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
if (fem-&gt;ops-&gt;integratejacobian) {ierr = (*fem-&gt;ops-&gt;integratejacobian)(fem, prob, jtype, fieldI, fieldJ, Ne, geom, coefficients, coefficients_t, probAux, coefficientsAux, elemMat);<A HREF="../Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</A>(ierr);}
<A HREF="../Sys/PetscFunctionReturn.html#PetscFunctionReturn">PetscFunctionReturn</A>(0);
}
<P>
#undef __FUNCT__
#define __FUNCT__ "PetscFEIntegrateBdJacobian"
/*C
PetscFEIntegrateBdJacobian - Produce the boundary element Jacobian for a chunk of elements by quadrature integration
<P>
Not collective
<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>fem          = The <A HREF="../DM/PetscFE.html#PetscFE">PetscFE</A> object for the field being integrated</B></TD><TD>- . prob         - The <A HREF="../DM/PetscDS.html#PetscDS">PetscDS</A> specifying the discretizations and continuum functions
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>fieldI       </B></TD><TD>- The test field being integrated
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>fieldJ       </B></TD><TD>- The basis field being integrated
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>Ne           </B></TD><TD>- The number of elements in the chunk
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>geom         </B></TD><TD>- The cell geometry for each cell in the chunk
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>coefficients </B></TD><TD>- The array of FEM basis coefficients for the elements for the Jacobian evaluation point
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>coefficients_t </B></TD><TD>- The array of FEM basis time derivative coefficients for the elements
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>probAux      </B></TD><TD>- The <A HREF="../DM/PetscDS.html#PetscDS">PetscDS</A> specifying the auxiliary discretizations
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>coefficientsAux </B></TD><TD>- The array of FEM auxiliary basis coefficients for the elements
</TD></TR></TABLE>
<P>
Output Parameter
<DT><B>elemMat              </B> -the element matrices for the Jacobian from each element
<br>
<P>
<H3><FONT COLOR="#CC3333">Note</FONT></H3>
<pre>
Loop over batch of elements (e):
</pre>
<pre>
  Loop over element matrix entries (f,fc,g,gc --&gt; i,j):
</pre>
<pre>
    Loop over quadrature points (q):
</pre>
<pre>
      Make u_q and gradU_q (loops over fields,Nb,Ncomp)
</pre>
<pre>
        elemMat[i,j] += \psi^{fc}_f(q) g0_{fc,gc}(u, \nabla u) \phi^{gc}_g(q)
</pre>
<pre>
                     + \psi^{fc}_f(q) \cdot g1_{fc,gc,dg}(u, \nabla u) \nabla\phi^{gc}_g(q)
</pre>
<pre>
                     + \nabla\psi^{fc}_f(q) \cdot g2_{fc,gc,df}(u, \nabla u) \phi^{gc}_g(q)
</pre>
<pre>
                     + \nabla\psi^{fc}_f(q) \cdot g3_{fc,gc,df,dg}(u, \nabla u) \nabla\phi^{gc}_g(q)
</pre>
*/
<A HREF="../Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</A> PetscFEIntegrateBdJacobian(<A HREF="../DM/PetscFE.html#PetscFE">PetscFE</A> fem, <A HREF="../DM/PetscDS.html#PetscDS">PetscDS</A> prob, <A HREF="../Sys/PetscInt.html#PetscInt">PetscInt</A> fieldI, <A HREF="../Sys/PetscInt.html#PetscInt">PetscInt</A> fieldJ, <A HREF="../Sys/PetscInt.html#PetscInt">PetscInt</A> Ne, PetscFECellGeom *geom,
const <A HREF="../Sys/PetscScalar.html#PetscScalar">PetscScalar</A> coefficients[], const <A HREF="../Sys/PetscScalar.html#PetscScalar">PetscScalar</A> coefficients_t[], <A HREF="../DM/PetscDS.html#PetscDS">PetscDS</A> probAux, const <A HREF="../Sys/PetscScalar.html#PetscScalar">PetscScalar</A> coefficientsAux[], <A HREF="../Sys/PetscScalar.html#PetscScalar">PetscScalar</A> elemMat[])
{
<A HREF="../Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</A> ierr;
<P>
<A HREF="../Sys/PetscFunctionBegin.html#PetscFunctionBegin">PetscFunctionBegin</A>;
PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
if (fem-&gt;ops-&gt;integratebdjacobian) {ierr = (*fem-&gt;ops-&gt;integratebdjacobian)(fem, prob, fieldI, fieldJ, Ne, geom, coefficients, coefficients_t, probAux, coefficientsAux, elemMat);<A HREF="../Sys/CHKERRQ.html#CHKERRQ">CHKERRQ</A>(ierr);}
<A HREF="../Sys/PetscFunctionReturn.html#PetscFunctionReturn">PetscFunctionReturn</A>(0);
}
<P>
#undef __FUNCT__
#define __FUNCT__ "PetscFERefine"
/*@
PetscFERefine - Create a "refined" <A HREF="../DM/PetscFE.html#PetscFE">PetscFE</A> object that refines the reference cell into smaller copies. This is typically used
to precondition a higher order method with a lower order method on a refined mesh having the same number of dofs (but more
sparsity). It is also used to create an interpolation between regularly refined meshes.
<P>
<H3><FONT COLOR="#CC3333">Input Parameter</FONT></H3>
<DT><B>fe </B> -The initial <A HREF="../DM/PetscFE.html#PetscFE">PetscFE</A>
<br>
<P>
<H3><FONT COLOR="#CC3333">Output Parameter</FONT></H3>
<DT><B>feRef </B> -The refined <A HREF="../DM/PetscFE.html#PetscFE">PetscFE</A>
<br>
<P>

<P>
<H3><FONT COLOR="#CC3333">See Also</FONT></H3>
 <A HREF="../DM/PetscFEType.html#PetscFEType">PetscFEType</A>, <A HREF="../DM/PetscFECreate.html#PetscFECreate">PetscFECreate</A>(), <A HREF="../DM/PetscFESetType.html#PetscFESetType">PetscFESetType</A>()
<BR><P><B><P><B><FONT COLOR="#CC3333">Level:</FONT></B>developer
<BR><FONT COLOR="#CC3333">Location:</FONT></B><A HREF="../../../src/dm/dt/interface/dtfe.c.html#PetscFEIntegrateJacobian">src/dm/dt/interface/dtfe.c</A>
<BR><A HREF="./index.html">Index of all DM 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>