xref: /petsc/src/snes/tutorials/network/water/water.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown static char help[] = "This example demonstrates the use of DMNetwork interface for solving a steady-state water network model.\n\
2*c4762a1bSJed Brown                       The water network equations follow those used for the package EPANET\n\
3*c4762a1bSJed Brown                       The data file format used is from the EPANET package (https://www.epa.gov/water-research/epanet).\n\
4*c4762a1bSJed Brown                       Run this program: mpiexec -n <n> ./water\n\\n";
5*c4762a1bSJed Brown 
6*c4762a1bSJed Brown /* T
7*c4762a1bSJed Brown    Concepts: DMNetwork
8*c4762a1bSJed Brown    Concepts: PETSc SNES solver
9*c4762a1bSJed Brown */
10*c4762a1bSJed Brown 
11*c4762a1bSJed Brown #include "water.h"
12*c4762a1bSJed Brown #include <petscdmnetwork.h>
13*c4762a1bSJed Brown 
14*c4762a1bSJed Brown int main(int argc,char ** argv)
15*c4762a1bSJed Brown {
16*c4762a1bSJed Brown   PetscErrorCode   ierr;
17*c4762a1bSJed Brown   char             waterdata_file[PETSC_MAX_PATH_LEN]="sample1.inp";
18*c4762a1bSJed Brown   WATERDATA        *waterdata;
19*c4762a1bSJed Brown   AppCtx_Water     appctx;
20*c4762a1bSJed Brown   PetscLogStage    stage1,stage2;
21*c4762a1bSJed Brown   PetscMPIInt      crank;
22*c4762a1bSJed Brown   DM               networkdm;
23*c4762a1bSJed Brown   PetscInt         *edgelist = NULL;
24*c4762a1bSJed Brown   PetscInt         nv,ne,i;
25*c4762a1bSJed Brown   const PetscInt   *vtx,*edges;
26*c4762a1bSJed Brown   Vec              X,F;
27*c4762a1bSJed Brown   SNES             snes;
28*c4762a1bSJed Brown   SNESConvergedReason reason;
29*c4762a1bSJed Brown 
30*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,"wateroptions",help);if (ierr) return ierr;
31*c4762a1bSJed Brown   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&crank);CHKERRQ(ierr);
32*c4762a1bSJed Brown 
33*c4762a1bSJed Brown   /* Create an empty network object */
34*c4762a1bSJed Brown   ierr = DMNetworkCreate(PETSC_COMM_WORLD,&networkdm);CHKERRQ(ierr);
35*c4762a1bSJed Brown 
36*c4762a1bSJed Brown   /* Register the components in the network */
37*c4762a1bSJed Brown   ierr = DMNetworkRegisterComponent(networkdm,"edgestruct",sizeof(struct _p_EDGE_Water),&appctx.compkey_edge);CHKERRQ(ierr);
38*c4762a1bSJed Brown   ierr = DMNetworkRegisterComponent(networkdm,"busstruct",sizeof(struct _p_VERTEX_Water),&appctx.compkey_vtx);CHKERRQ(ierr);
39*c4762a1bSJed Brown 
40*c4762a1bSJed Brown   ierr = PetscLogStageRegister("Read Data",&stage1);CHKERRQ(ierr);
41*c4762a1bSJed Brown   PetscLogStagePush(stage1);
42*c4762a1bSJed Brown   ierr = PetscNew(&waterdata);CHKERRQ(ierr);
43*c4762a1bSJed Brown 
44*c4762a1bSJed Brown   /* READ THE DATA */
45*c4762a1bSJed Brown   if (!crank) {
46*c4762a1bSJed Brown     /* READ DATA. Only rank 0 reads the data */
47*c4762a1bSJed Brown     ierr = PetscOptionsGetString(NULL,NULL,"-waterdata",waterdata_file,PETSC_MAX_PATH_LEN-1,NULL);CHKERRQ(ierr);
48*c4762a1bSJed Brown     ierr = WaterReadData(waterdata,waterdata_file);CHKERRQ(ierr);
49*c4762a1bSJed Brown 
50*c4762a1bSJed Brown     ierr = PetscCalloc1(2*waterdata->nedge,&edgelist);CHKERRQ(ierr);
51*c4762a1bSJed Brown     ierr = GetListofEdges_Water(waterdata,edgelist);CHKERRQ(ierr);
52*c4762a1bSJed Brown   }
53*c4762a1bSJed Brown   PetscLogStagePop();
54*c4762a1bSJed Brown 
55*c4762a1bSJed Brown   ierr = PetscLogStageRegister("Create network",&stage2);CHKERRQ(ierr);
56*c4762a1bSJed Brown   PetscLogStagePush(stage2);
57*c4762a1bSJed Brown 
58*c4762a1bSJed Brown   /* Set numbers of nodes and edges */
59*c4762a1bSJed Brown   ierr = DMNetworkSetSizes(networkdm,1,&waterdata->nvertex,&waterdata->nedge,0,NULL);CHKERRQ(ierr);
60*c4762a1bSJed Brown   if (!crank) {
61*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_SELF,"water nvertices %D, nedges %D\n",waterdata->nvertex,waterdata->nedge);CHKERRQ(ierr);
62*c4762a1bSJed Brown   }
63*c4762a1bSJed Brown 
64*c4762a1bSJed Brown   /* Add edge connectivity */
65*c4762a1bSJed Brown   ierr = DMNetworkSetEdgeList(networkdm,&edgelist,NULL);CHKERRQ(ierr);
66*c4762a1bSJed Brown 
67*c4762a1bSJed Brown   /* Set up the network layout */
68*c4762a1bSJed Brown   ierr = DMNetworkLayoutSetUp(networkdm);CHKERRQ(ierr);
69*c4762a1bSJed Brown 
70*c4762a1bSJed Brown   if (!crank) {
71*c4762a1bSJed Brown     ierr = PetscFree(edgelist);CHKERRQ(ierr);
72*c4762a1bSJed Brown   }
73*c4762a1bSJed Brown 
74*c4762a1bSJed Brown   /* ADD VARIABLES AND COMPONENTS FOR THE NETWORK */
75*c4762a1bSJed Brown   ierr = DMNetworkGetSubnetworkInfo(networkdm,0,&nv,&ne,&vtx,&edges);CHKERRQ(ierr);
76*c4762a1bSJed Brown 
77*c4762a1bSJed Brown   for (i = 0; i < ne; i++) {
78*c4762a1bSJed Brown     ierr = DMNetworkAddComponent(networkdm,edges[i],appctx.compkey_edge,&waterdata->edge[i]);CHKERRQ(ierr);
79*c4762a1bSJed Brown   }
80*c4762a1bSJed Brown 
81*c4762a1bSJed Brown   for (i = 0; i < nv; i++) {
82*c4762a1bSJed Brown     ierr = DMNetworkAddComponent(networkdm,vtx[i],appctx.compkey_vtx,&waterdata->vertex[i]);CHKERRQ(ierr);
83*c4762a1bSJed Brown     /* Add number of variables */
84*c4762a1bSJed Brown     ierr = DMNetworkAddNumVariables(networkdm,vtx[i],1);CHKERRQ(ierr);
85*c4762a1bSJed Brown   }
86*c4762a1bSJed Brown 
87*c4762a1bSJed Brown   /* Set up DM for use */
88*c4762a1bSJed Brown   ierr = DMSetUp(networkdm);CHKERRQ(ierr);
89*c4762a1bSJed Brown 
90*c4762a1bSJed Brown   if (!crank) {
91*c4762a1bSJed Brown     ierr = PetscFree(waterdata->vertex);CHKERRQ(ierr);
92*c4762a1bSJed Brown     ierr = PetscFree(waterdata->edge);CHKERRQ(ierr);
93*c4762a1bSJed Brown   }
94*c4762a1bSJed Brown   ierr = PetscFree(waterdata);CHKERRQ(ierr);
95*c4762a1bSJed Brown 
96*c4762a1bSJed Brown   /* Distribute networkdm to multiple processes */
97*c4762a1bSJed Brown   ierr = DMNetworkDistribute(&networkdm,0);CHKERRQ(ierr);
98*c4762a1bSJed Brown 
99*c4762a1bSJed Brown   PetscLogStagePop();
100*c4762a1bSJed Brown 
101*c4762a1bSJed Brown   ierr = DMCreateGlobalVector(networkdm,&X);CHKERRQ(ierr);
102*c4762a1bSJed Brown   ierr = VecDuplicate(X,&F);CHKERRQ(ierr);
103*c4762a1bSJed Brown 
104*c4762a1bSJed Brown   /* HOOK UP SOLVER */
105*c4762a1bSJed Brown   ierr = SNESCreate(PETSC_COMM_WORLD,&snes);CHKERRQ(ierr);
106*c4762a1bSJed Brown   ierr = SNESSetDM(snes,networkdm);CHKERRQ(ierr);
107*c4762a1bSJed Brown   ierr = SNESSetOptionsPrefix(snes,"water_");CHKERRQ(ierr);
108*c4762a1bSJed Brown   ierr = SNESSetFunction(snes,F,WaterFormFunction,NULL);CHKERRQ(ierr);
109*c4762a1bSJed Brown   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
110*c4762a1bSJed Brown 
111*c4762a1bSJed Brown   ierr = WaterSetInitialGuess(networkdm,X);CHKERRQ(ierr);
112*c4762a1bSJed Brown   /* ierr = VecView(X,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */
113*c4762a1bSJed Brown 
114*c4762a1bSJed Brown   ierr = SNESSolve(snes,NULL,X);CHKERRQ(ierr);
115*c4762a1bSJed Brown   ierr = SNESGetConvergedReason(snes,&reason);CHKERRQ(ierr);
116*c4762a1bSJed Brown   if (reason < 0) {
117*c4762a1bSJed Brown     SETERRQ(PETSC_COMM_SELF,0,"No solution found for the water network");
118*c4762a1bSJed Brown   }
119*c4762a1bSJed Brown   /* ierr = VecView(X,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */
120*c4762a1bSJed Brown 
121*c4762a1bSJed Brown   ierr = VecDestroy(&X);CHKERRQ(ierr);
122*c4762a1bSJed Brown   ierr = VecDestroy(&F);CHKERRQ(ierr);
123*c4762a1bSJed Brown   ierr = SNESDestroy(&snes);CHKERRQ(ierr);
124*c4762a1bSJed Brown   ierr = DMDestroy(&networkdm);CHKERRQ(ierr);
125*c4762a1bSJed Brown   ierr = PetscFinalize();
126*c4762a1bSJed Brown   return ierr;
127*c4762a1bSJed Brown }
128*c4762a1bSJed Brown 
129*c4762a1bSJed Brown /*TEST
130*c4762a1bSJed Brown 
131*c4762a1bSJed Brown    build:
132*c4762a1bSJed Brown       depends: waterreaddata.c waterfunctions.c
133*c4762a1bSJed Brown       requires: !complex double define(PETSC_HAVE_ATTRIBUTEALIGNED)
134*c4762a1bSJed Brown 
135*c4762a1bSJed Brown    test:
136*c4762a1bSJed Brown       args: -water_snes_converged_reason -options_left no
137*c4762a1bSJed Brown       localrunfiles: wateroptions sample1.inp
138*c4762a1bSJed Brown       output_file: output/water.out
139*c4762a1bSJed Brown       requires: double !complex define(PETSC_HAVE_ATTRIBUTEALIGNED)
140*c4762a1bSJed Brown 
141*c4762a1bSJed Brown    test:
142*c4762a1bSJed Brown       suffix: 2
143*c4762a1bSJed Brown       nsize: 3
144*c4762a1bSJed Brown       args: -water_snes_converged_reason -options_left no
145*c4762a1bSJed Brown       localrunfiles: wateroptions sample1.inp
146*c4762a1bSJed Brown       output_file: output/water.out
147*c4762a1bSJed Brown       requires: double !complex define(PETSC_HAVE_ATTRIBUTEALIGNED)
148*c4762a1bSJed Brown 
149*c4762a1bSJed Brown TEST*/
150