xref: /petsc/src/ts/tutorials/network/wash.h (revision df4cd43f92eaa320656440c40edb1046daee8f75)
1 #ifndef WASH_H
2 #define WASH_H
3 
4 #include <petscdmnetwork.h>
5 #include "pipe.h"
6 
7 typedef enum {
8   NONE,
9   JUNCTION  = 1,
10   RESERVOIR = 2,
11   VALVE     = 3,
12   DEMAND    = 4,
13   INFLOW    = 5,
14   STAGE     = 6,
15   TANK      = 7
16 } VertexType;
17 
18 typedef struct {
19   PetscInt    rid;  /*reservoir id*/
20   PetscScalar hres; /*Reservoir water column*/
21 } Reservoir;
22 
23 typedef struct {
24   PetscInt    vid; /*valve id*/
25   PetscScalar tau; /*valve aperture*/
26   PetscScalar cdag;
27   PetscScalar qout;
28 } Valve;
29 
30 /* junction              */
31 /*-----------------------*/
32 struct _p_Junction {
33   PetscInt   id;  /* global index */
34   PetscInt   tag; /* external id */
35   VertexType type;
36   PetscInt   isEnd;                 /* -1: left end; 0: not an end; 1: right end */
37   PetscInt   nedges_in, nedges_out; /* number of connected in/out edges */
38   Mat       *jacobian;
39   PetscReal  latitude, longitude; /* GPS data */
40 
41   /* boundary data structures */
42   Reservoir reservoir;
43   Valve     valve;
44 } PETSC_ATTRIBUTEALIGNED(PetscMax(sizeof(double), sizeof(PetscScalar)));
45 typedef struct _p_Junction *Junction;
46 
47 extern PetscErrorCode JunctionCreateJacobian(DM, PetscInt, Mat *, Mat *[]);
48 extern PetscErrorCode JunctionDestroyJacobian(DM, PetscInt, Junction);
49 
50 /* wash                   */
51 /*------------------------*/
52 struct _p_Wash {
53   MPI_Comm  comm;
54   PetscInt  nedge, nvertex;    /* local number of components */
55   PetscInt  Nedge, Nvertex;    /* global number of components */
56   PetscInt *edgelist;          /* local edge list */
57   Vec       localX, localXdot; /* vectors used in local function evaluation */
58   PetscInt  nnodes_loc;        /* num of global and local nodes */
59 
60   /* Junction */
61   Junction  junction;
62   PetscInt *vtype;
63 
64   /* Pipe */
65   Pipe        pipe;
66   PetscScalar Q0, H0, QL, HL; /* left and right boundary conditions for wash-network (not individual pipe) */
67 
68   /* Events */
69   PetscInt close_valve;
70 } PETSC_ATTRIBUTEALIGNED(PetscMax(sizeof(double), sizeof(PetscScalar)));
71 typedef struct _p_Wash *Wash;
72 
73 extern PetscErrorCode WashNetworkCreate(MPI_Comm, PetscInt, Wash *);
74 extern PetscErrorCode WashNetworkCleanUp(Wash);
75 #endif
76