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
|
<center><a href="https://gitlab.com/petsc/petsc/-/blob/966382dc56242773704ef5f5cee7aa2db3ebc577/include/petscdstypes.h">Actual source code: petscdstypes.h</a></center><br>
<html>
<head>
<title></title>
<meta name="generator" content="c2html 0.9.6">
<meta name="date" content="2025-04-30T18:14:50+00:00">
</head>
<body bgcolor="#FFFFFF">
<pre width=80>
<a name="line1"> 1: </a><font color="#A020F0">#pragma once</font>
<a name="line3"> 3: </a>#include <A href="../include/petscdmlabel.h.html"><petscdmlabel.h></A>
<a name="line5"> 5: </a><font color="#B22222">/* MANSEC = <a href="../manualpages/DM/DM.html">DM</a> */</font>
<a name="line6"> 6: </a><font color="#B22222">/* SUBMANSEC = DT */</font>
<a name="line8"> 8: </a><font color="#B22222">/*S</font>
<a name="line9"> 9: </a><font color="#B22222"> <a href="../manualpages/DT/PetscDS.html">PetscDS</a> - PETSc object that manages a discrete system, which is a set of discretizations + continuum equations from a `<a href="../manualpages/DT/PetscWeakForm.html">PetscWeakForm</a>`</font>
<a name="line11"> 11: </a><font color="#B22222"> Level: intermediate</font>
<a name="line13"> 13: </a><font color="#B22222">.seealso: `<a href="../manualpages/DT/PetscDSCreate.html">PetscDSCreate</a>()`, `<a href="../manualpages/DT/PetscDSSetType.html">PetscDSSetType</a>()`, `<a href="../manualpages/DT/PetscDSType.html">PetscDSType</a>`, `<a href="../manualpages/DT/PetscWeakForm.html">PetscWeakForm</a>`, `<a href="../manualpages/FE/PetscFECreate.html">PetscFECreate</a>()`, `<a href="../manualpages/FV/PetscFVCreate.html">PetscFVCreate</a>()`</font>
<a name="line14"> 14: </a><font color="#B22222">S*/</font>
<a name="line15"> 15: </a><font color="#4169E1">typedef struct _p_PetscDS *<a href="../manualpages/DT/PetscDS.html">PetscDS</a>;</font>
<a name="line17"> 17: </a><font color="#B22222">/*S</font>
<a name="line18"> 18: </a><font color="#B22222"> <a href="../manualpages/DT/PetscWeakForm.html">PetscWeakForm</a> - PETSc object that manages a sets of pointwise functions defining a system of equations</font>
<a name="line20"> 20: </a><font color="#B22222"> Level: intermediate</font>
<a name="line22"> 22: </a><font color="#B22222">.seealso: `<a href="../manualpages/DT/PetscWeakFormCreate.html">PetscWeakFormCreate</a>()`, `<a href="../manualpages/DT/PetscDS.html">PetscDS</a>`, `<a href="../manualpages/FE/PetscFECreate.html">PetscFECreate</a>()`, `<a href="../manualpages/FV/PetscFVCreate.html">PetscFVCreate</a>()`</font>
<a name="line23"> 23: </a><font color="#B22222">S*/</font>
<a name="line24"> 24: </a><font color="#4169E1">typedef struct _p_PetscWeakForm *<a href="../manualpages/DT/PetscWeakForm.html">PetscWeakForm</a>;</font>
<a name="line26"> 26: </a><font color="#B22222">/*S</font>
<a name="line27"> 27: </a><font color="#B22222"> <a href="../manualpages/DT/PetscFormKey.html">PetscFormKey</a> - This key indicates how to use a set of pointwise functions defining part of a system of equations</font>
<a name="line29"> 29: </a><font color="#B22222"> The subdomain on which to integrate is specified by (label, value), the test function field by (field), and the</font>
<a name="line30"> 30: </a><font color="#B22222"> piece of the equation by (part). For example, LHS = 0 and RHS = 1 in IMEX methods. More pieces can be present for</font>
<a name="line31"> 31: </a><font color="#B22222"> operator splitting methods.</font>
<a name="line33"> 33: </a><font color="#B22222"> Level: intermediate</font>
<a name="line35"> 35: </a><font color="#B22222"> Note:</font>
<a name="line36"> 36: </a><font color="#B22222"> This is a struct, not a `<a href="../manualpages/Sys/PetscObject.html">PetscObject</a>`</font>
<a name="line38"> 38: </a><font color="#B22222">.seealso: `<a href="../manualpages/SNES/DMPlexSNESComputeResidualFEM.html">DMPlexSNESComputeResidualFEM</a>()`, `<a href="../manualpages/SNES/DMPlexSNESComputeJacobianFEM.html">DMPlexSNESComputeJacobianFEM</a>()`, `<a href="../manualpages/SNES/DMPlexSNESComputeBoundaryFEM.html">DMPlexSNESComputeBoundaryFEM</a>()`</font>
<a name="line39"> 39: </a><font color="#B22222">S*/</font>
<a name="line40"> 40: </a><font color="#4169E1">typedef</font> <font color="#4169E1">struct</font> {
<a name="line41"> 41: </a> <a href="../manualpages/DM/DMLabel.html">DMLabel</a> label; <font color="#B22222">/* The (label, value) select a subdomain */</font>
<a name="line42"> 42: </a> <a href="../manualpages/Sys/PetscInt.html">PetscInt</a> value;
<a name="line43"> 43: </a> <a href="../manualpages/Sys/PetscInt.html">PetscInt</a> field; <font color="#B22222">/* Selects the field for the test function */</font>
<a name="line44"> 44: </a> <a href="../manualpages/Sys/PetscInt.html">PetscInt</a> part; <font color="#B22222">/* Selects the equation part. For example, LHS = 0 and RHS = 1 in IMEX methods. More pieces can be present for operator splitting methods. */</font>
<a name="line45"> 45: </a>} <a href="../manualpages/DT/PetscFormKey.html">PetscFormKey</a>;
<a name="line47"> 47: </a><font color="#B22222">/*E</font>
<a name="line48"> 48: </a><font color="#B22222"> <a href="../manualpages/DT/PetscWeakFormKind.html">PetscWeakFormKind</a> - The kind of weak form. The specific forms are given in the documentation for the integraton functions.</font>
<a name="line50"> 50: </a><font color="#B22222"> Values:</font>
<a name="line51"> 51: </a><font color="#B22222">+ OBJECTIVE - Objective form</font>
<a name="line52"> 52: </a><font color="#B22222">. F0, F1 - Residual forms</font>
<a name="line53"> 53: </a><font color="#B22222">. G0, G1, G2, G3 - Jacobian forms</font>
<a name="line54"> 54: </a><font color="#B22222">. GP0, GP1, GP2, GP3 - Jacobian preconditioner matrix forms</font>
<a name="line55"> 55: </a><font color="#B22222">. GT0, GT1, GT2, GT3 - Dynamic Jacobian matrix forms</font>
<a name="line56"> 56: </a><font color="#B22222">. BDF0, BDF1 - Boundary Residual forms</font>
<a name="line57"> 57: </a><font color="#B22222">. BDG0, BDG1, BDG2, BDG3 - Jacobian forms</font>
<a name="line58"> 58: </a><font color="#B22222">. BDGP0, BDGP1, BDGP2, BDGP3 - Jacobian preconditioner matrix forms</font>
<a name="line59"> 59: </a><font color="#B22222">. R - Riemann solver</font>
<a name="line60"> 60: </a><font color="#B22222">- CEED - libCEED QFunction</font>
<a name="line62"> 62: </a><font color="#B22222"> Level: beginner</font>
<a name="line64"> 64: </a><font color="#B22222">.seealso: `<a href="../manualpages/DT/PetscWeakForm.html">PetscWeakForm</a>`, `<a href="../manualpages/FE/PetscFEIntegrateResidual.html">PetscFEIntegrateResidual</a>()`, `<a href="../manualpages/FE/PetscFEIntegrateJacobian.html">PetscFEIntegrateJacobian</a>()`, `<a href="../manualpages/FE/PetscFEIntegrateBdResidual.html">PetscFEIntegrateBdResidual</a>()`, `<a href="../manualpages/FE/PetscFEIntegrateBdJacobian.html">PetscFEIntegrateBdJacobian</a>()`,</font>
<a name="line65"> 65: </a><font color="#B22222"> `<a href="../manualpages/FV/PetscFVIntegrateRHSFunction.html">PetscFVIntegrateRHSFunction</a>()`, `PetscWeakFormSetIndexResidual()`, `PetscWeakFormClearIndex()`</font>
<a name="line66"> 66: </a><font color="#B22222">E*/</font>
<a name="line67"> 67: </a><font color="#4169E1">typedef</font> <font color="#4169E1">enum</font> {
<a name="line68"> 68: </a> <a href="../manualpages/DT/PetscWeakFormKind.html">PETSC_WF_OBJECTIVE</a>,
<a name="line69"> 69: </a> <a href="../manualpages/DT/PetscWeakFormKind.html">PETSC_WF_F0</a>,
<a name="line70"> 70: </a> <a href="../manualpages/DT/PetscWeakFormKind.html">PETSC_WF_F1</a>,
<a name="line71"> 71: </a> <a href="../manualpages/DT/PetscWeakFormKind.html">PETSC_WF_G0</a>,
<a name="line72"> 72: </a> <a href="../manualpages/DT/PetscWeakFormKind.html">PETSC_WF_G1</a>,
<a name="line73"> 73: </a> <a href="../manualpages/DT/PetscWeakFormKind.html">PETSC_WF_G2</a>,
<a name="line74"> 74: </a> <a href="../manualpages/DT/PetscWeakFormKind.html">PETSC_WF_G3</a>,
<a name="line75"> 75: </a> <a href="../manualpages/DT/PetscWeakFormKind.html">PETSC_WF_GP0</a>,
<a name="line76"> 76: </a> <a href="../manualpages/DT/PetscWeakFormKind.html">PETSC_WF_GP1</a>,
<a name="line77"> 77: </a> <a href="../manualpages/DT/PetscWeakFormKind.html">PETSC_WF_GP2</a>,
<a name="line78"> 78: </a> <a href="../manualpages/DT/PetscWeakFormKind.html">PETSC_WF_GP3</a>,
<a name="line79"> 79: </a> <a href="../manualpages/DT/PetscWeakFormKind.html">PETSC_WF_GT0</a>,
<a name="line80"> 80: </a> <a href="../manualpages/DT/PetscWeakFormKind.html">PETSC_WF_GT1</a>,
<a name="line81"> 81: </a> <a href="../manualpages/DT/PetscWeakFormKind.html">PETSC_WF_GT2</a>,
<a name="line82"> 82: </a> <a href="../manualpages/DT/PetscWeakFormKind.html">PETSC_WF_GT3</a>,
<a name="line83"> 83: </a> <a href="../manualpages/DT/PetscWeakFormKind.html">PETSC_WF_BDF0</a>,
<a name="line84"> 84: </a> <a href="../manualpages/DT/PetscWeakFormKind.html">PETSC_WF_BDF1</a>,
<a name="line85"> 85: </a> <a href="../manualpages/DT/PetscWeakFormKind.html">PETSC_WF_BDG0</a>,
<a name="line86"> 86: </a> <a href="../manualpages/DT/PetscWeakFormKind.html">PETSC_WF_BDG1</a>,
<a name="line87"> 87: </a> <a href="../manualpages/DT/PetscWeakFormKind.html">PETSC_WF_BDG2</a>,
<a name="line88"> 88: </a> <a href="../manualpages/DT/PetscWeakFormKind.html">PETSC_WF_BDG3</a>,
<a name="line89"> 89: </a> <a href="../manualpages/DT/PetscWeakFormKind.html">PETSC_WF_BDGP0</a>,
<a name="line90"> 90: </a> <a href="../manualpages/DT/PetscWeakFormKind.html">PETSC_WF_BDGP1</a>,
<a name="line91"> 91: </a> <a href="../manualpages/DT/PetscWeakFormKind.html">PETSC_WF_BDGP2</a>,
<a name="line92"> 92: </a> <a href="../manualpages/DT/PetscWeakFormKind.html">PETSC_WF_BDGP3</a>,
<a name="line93"> 93: </a> <a href="../manualpages/DT/PetscWeakFormKind.html">PETSC_WF_R</a>,
<a name="line94"> 94: </a> <a href="../manualpages/DT/PetscWeakFormKind.html">PETSC_WF_CEED</a>,
<a name="line95"> 95: </a> <a href="../manualpages/DT/PetscWeakFormKind.html">PETSC_NUM_WF</a>
<a name="line96"> 96: </a>} <a href="../manualpages/DT/PetscWeakFormKind.html">PetscWeakFormKind</a>;
<a name="line97"> 97: </a>PETSC_EXTERN const char *const PetscWeakFormKinds[];
</pre>
</body>
</html>
|