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