xref: /petsc/src/ts/tutorials/network/pipe.h (revision 9dd11ecf0918283bb567d8b33a92f53ac4ea7840)
1 #pragma once
2 
3 #define GRAV                9.806
4 #define PIPE_CHARACTERISTIC 10000000.0
5 
6 #include <petsc.h>
7 
8 typedef struct {
9   PetscScalar q; /* flow rate */
10   PetscScalar h; /* pressure */
11 } PipeField;
12 
13 typedef struct {
14   PetscScalar Q0, H0; /* boundary values in upstream */
15   PetscScalar QL, HL; /* boundary values in downstream */
16 } PipeBoundary;
17 
18 /* pipe                 */
19 /*----------------------*/
20 struct _p_Pipe {
21   /* identification variables */
22   PetscInt id;
23   PetscInt networkid; /* which network this pipe belongs */
24 
25   /* solver objects */
26   Vec        x;
27   PipeField *xold;
28   PetscReal  dt;
29   DM         da;
30   PetscInt   nnodes; /* number of nodes in da discretization */
31   Mat       *jacobian;
32 
33   /* physics */
34   PetscReal    length; /* pipe length */
35   PetscReal    a;      /* natural flow speed */
36   PetscReal    fric;   /* friction */
37   PetscReal    D;      /* diameter */
38   PetscReal    A;      /* area of cross section */
39   PetscReal    R;
40   PetscReal    rad;
41   PipeBoundary boundary; /* boundary conditions for H and Q */
42 } PETSC_ATTRIBUTEALIGNED(PetscMax(sizeof(double), sizeof(PetscScalar)));
43 
44 typedef struct _p_Pipe *Pipe;
45 
46 extern PetscErrorCode PipeCreate(MPI_Comm, Pipe *);
47 extern PetscErrorCode PipeDestroy(Pipe *);
48 extern PetscErrorCode PipeSetParameters(Pipe, PetscReal, PetscReal, PetscReal, PetscReal);
49 extern PetscErrorCode PipeSetUp(Pipe);
50 extern PetscErrorCode PipeCreateJacobian(Pipe, Mat *, Mat *[]);
51 extern PetscErrorCode PipeDestroyJacobian(Pipe);
52 
53 extern PetscErrorCode PipeComputeSteadyState(Pipe, PetscScalar, PetscScalar);
54 extern PetscErrorCode PipeIFunctionLocal(DMDALocalInfo *, PetscReal, PipeField *, PipeField *, PipeField *, Pipe);
55 extern PetscErrorCode PipeIFunctionLocal_Lax(DMDALocalInfo *, PetscReal, PipeField *, PipeField *, PetscScalar *, Pipe);
56 extern PetscErrorCode PipeRHSFunctionLocal(DMDALocalInfo *, PetscReal, PipeField *, PetscScalar *, Pipe);
57 extern PetscErrorCode PipeMonitor(TS, PetscInt, PetscReal, Vec, void *);
58 
59 extern PetscErrorCode PipeCreateJacobian(Pipe, Mat *, Mat *[]);
60 extern PetscErrorCode PipeDestroyJacobian(Pipe);
61