File: PetscFEIntegrateJacobian.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 (89 lines) | stat: -rw-r--r-- 5,788 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
<!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.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/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" 
<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, PetscFEGeom *cgeom,
                                        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/PetscReal.html#PetscReal">PetscReal</A> t, <A HREF="../Sys/PetscReal.html#PetscReal">PetscReal</A> u_tshift, <A HREF="../Sys/PetscScalar.html#PetscScalar">PetscScalar</A> elemMat[])
</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>cgeom        </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>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>t            </B></TD><TD>- The time
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>u_tShift     </B></TD><TD>- A multiplier for the dF/du_t term (as opposed to the dF/du term)
</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>

<P>
<H3><FONT COLOR="#CC3333">See Also</FONT></H3>
 <A HREF="../DM/PetscFEIntegrateResidual.html#PetscFEIntegrateResidual">PetscFEIntegrateResidual</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/fe/interface/fe.c.html#PetscFEIntegrateJacobian">src/dm/dt/fe/interface/fe.c</A>
<P><H3><FONT COLOR="CC3333">Implementations</FONT></H3>src/dm/dt/fe/impls/basic/febasic.c:577:PetscErrorCode PetscFEIntegrateJacobian_Basic(PetscFE fem, PetscDS prob, PetscFEJacobianType jtype, PetscInt fieldI, PetscInt fieldJ, PetscInt Ne, PetscFEGeom *geom,
<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>