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