1*c4762a1bSJed Brown static char help[] = "This example demonstrates the use of DMNetwork interface with subnetworks for solving a coupled nonlinear \n\ 2*c4762a1bSJed Brown electric power grid and water pipe problem.\n\ 3*c4762a1bSJed Brown The available solver options are in the ex1options file \n\ 4*c4762a1bSJed Brown and the data files are in the datafiles of subdirectories.\n\ 5*c4762a1bSJed Brown This example shows the use of subnetwork feature in DMNetwork. \n\ 6*c4762a1bSJed Brown Run this program: mpiexec -n <n> ./ex1 \n\\n"; 7*c4762a1bSJed Brown 8*c4762a1bSJed Brown /* T 9*c4762a1bSJed Brown Concepts: DMNetwork 10*c4762a1bSJed Brown Concepts: PETSc SNES solver 11*c4762a1bSJed Brown */ 12*c4762a1bSJed Brown 13*c4762a1bSJed Brown #include "power/power.h" 14*c4762a1bSJed Brown #include "water/water.h" 15*c4762a1bSJed Brown 16*c4762a1bSJed Brown typedef struct{ 17*c4762a1bSJed Brown UserCtx_Power appctx_power; 18*c4762a1bSJed Brown AppCtx_Water appctx_water; 19*c4762a1bSJed Brown PetscInt subsnes_id; /* snes solver id */ 20*c4762a1bSJed Brown PetscInt it; /* iteration number */ 21*c4762a1bSJed Brown Vec localXold; /* store previous solution, used by FormFunction_Dummy() */ 22*c4762a1bSJed Brown } UserCtx; 23*c4762a1bSJed Brown 24*c4762a1bSJed Brown PetscErrorCode UserMonitor(SNES snes,PetscInt its,PetscReal fnorm ,void *appctx) 25*c4762a1bSJed Brown { 26*c4762a1bSJed Brown PetscErrorCode ierr; 27*c4762a1bSJed Brown UserCtx *user = (UserCtx*)appctx; 28*c4762a1bSJed Brown Vec X,localXold=user->localXold; 29*c4762a1bSJed Brown DM networkdm; 30*c4762a1bSJed Brown PetscMPIInt rank; 31*c4762a1bSJed Brown MPI_Comm comm; 32*c4762a1bSJed Brown 33*c4762a1bSJed Brown PetscFunctionBegin; 34*c4762a1bSJed Brown ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr); 35*c4762a1bSJed Brown ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 36*c4762a1bSJed Brown #if 0 37*c4762a1bSJed Brown if (!rank) { 38*c4762a1bSJed Brown PetscInt subsnes_id=user->subsnes_id; 39*c4762a1bSJed Brown if (subsnes_id == 2) { 40*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF," it %d, subsnes_id %d, fnorm %g\n",user->it,user->subsnes_id,fnorm);CHKERRQ(ierr); 41*c4762a1bSJed Brown } else { 42*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF," subsnes_id %d, fnorm %g\n",user->subsnes_id,fnorm);CHKERRQ(ierr); 43*c4762a1bSJed Brown } 44*c4762a1bSJed Brown } 45*c4762a1bSJed Brown #endif 46*c4762a1bSJed Brown ierr = SNESGetSolution(snes,&X);CHKERRQ(ierr); 47*c4762a1bSJed Brown ierr = SNESGetDM(snes,&networkdm);CHKERRQ(ierr); 48*c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(networkdm,X,INSERT_VALUES,localXold);CHKERRQ(ierr); 49*c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(networkdm,X,INSERT_VALUES,localXold);CHKERRQ(ierr); 50*c4762a1bSJed Brown PetscFunctionReturn(0); 51*c4762a1bSJed Brown } 52*c4762a1bSJed Brown 53*c4762a1bSJed Brown PetscErrorCode FormJacobian_subPower(SNES snes,Vec X, Mat J,Mat Jpre,void *appctx) 54*c4762a1bSJed Brown { 55*c4762a1bSJed Brown PetscErrorCode ierr; 56*c4762a1bSJed Brown DM networkdm; 57*c4762a1bSJed Brown Vec localX; 58*c4762a1bSJed Brown PetscInt nv,ne,i,j,offset,nvar,row; 59*c4762a1bSJed Brown const PetscInt *vtx,*edges; 60*c4762a1bSJed Brown PetscBool ghostvtex; 61*c4762a1bSJed Brown PetscScalar one = 1.0; 62*c4762a1bSJed Brown PetscMPIInt rank; 63*c4762a1bSJed Brown MPI_Comm comm; 64*c4762a1bSJed Brown 65*c4762a1bSJed Brown PetscFunctionBegin; 66*c4762a1bSJed Brown ierr = SNESGetDM(snes,&networkdm);CHKERRQ(ierr); 67*c4762a1bSJed Brown ierr = DMGetLocalVector(networkdm,&localX);CHKERRQ(ierr); 68*c4762a1bSJed Brown 69*c4762a1bSJed Brown ierr = PetscObjectGetComm((PetscObject)networkdm,&comm);CHKERRQ(ierr); 70*c4762a1bSJed Brown ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 71*c4762a1bSJed Brown 72*c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(networkdm,X,INSERT_VALUES,localX);CHKERRQ(ierr); 73*c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(networkdm,X,INSERT_VALUES,localX);CHKERRQ(ierr); 74*c4762a1bSJed Brown 75*c4762a1bSJed Brown ierr = MatZeroEntries(J);CHKERRQ(ierr); 76*c4762a1bSJed Brown 77*c4762a1bSJed Brown /* Power subnetwork: copied from power/FormJacobian_Power() */ 78*c4762a1bSJed Brown ierr = DMNetworkGetSubnetworkInfo(networkdm,0,&nv,&ne,&vtx,&edges);CHKERRQ(ierr); 79*c4762a1bSJed Brown ierr = FormJacobian_Power_private(networkdm,localX,J,nv,ne,vtx,edges,appctx);CHKERRQ(ierr); 80*c4762a1bSJed Brown 81*c4762a1bSJed Brown /* Water subnetwork: Identity */ 82*c4762a1bSJed Brown ierr = DMNetworkGetSubnetworkInfo(networkdm,1,&nv,&ne,&vtx,&edges);CHKERRQ(ierr); 83*c4762a1bSJed Brown for (i=0; i<nv; i++) { 84*c4762a1bSJed Brown ierr = DMNetworkIsGhostVertex(networkdm,vtx[i],&ghostvtex);CHKERRQ(ierr); 85*c4762a1bSJed Brown if (ghostvtex) continue; 86*c4762a1bSJed Brown 87*c4762a1bSJed Brown ierr = DMNetworkGetVariableGlobalOffset(networkdm,vtx[i],&offset);CHKERRQ(ierr); 88*c4762a1bSJed Brown ierr = DMNetworkGetNumVariables(networkdm,vtx[i],&nvar);CHKERRQ(ierr); 89*c4762a1bSJed Brown for (j=0; j<nvar; j++) { 90*c4762a1bSJed Brown row = offset + j; 91*c4762a1bSJed Brown ierr = MatSetValues(J,1,&row,1,&row,&one,ADD_VALUES);CHKERRQ(ierr); 92*c4762a1bSJed Brown } 93*c4762a1bSJed Brown } 94*c4762a1bSJed Brown ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 95*c4762a1bSJed Brown ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 96*c4762a1bSJed Brown 97*c4762a1bSJed Brown ierr = DMRestoreLocalVector(networkdm,&localX);CHKERRQ(ierr); 98*c4762a1bSJed Brown PetscFunctionReturn(0); 99*c4762a1bSJed Brown } 100*c4762a1bSJed Brown 101*c4762a1bSJed Brown /* Dummy equation localF(X) = localX - localXold */ 102*c4762a1bSJed Brown PetscErrorCode FormFunction_Dummy(DM networkdm,Vec localX, Vec localF,PetscInt nv,PetscInt ne,const PetscInt* vtx,const PetscInt* edges,void* appctx) 103*c4762a1bSJed Brown { 104*c4762a1bSJed Brown PetscErrorCode ierr; 105*c4762a1bSJed Brown const PetscScalar *xarr,*xoldarr; 106*c4762a1bSJed Brown PetscScalar *farr; 107*c4762a1bSJed Brown PetscInt i,j,offset,nvar; 108*c4762a1bSJed Brown PetscBool ghostvtex; 109*c4762a1bSJed Brown UserCtx *user = (UserCtx*)appctx; 110*c4762a1bSJed Brown Vec localXold = user->localXold; 111*c4762a1bSJed Brown 112*c4762a1bSJed Brown PetscFunctionBegin; 113*c4762a1bSJed Brown ierr = VecGetArrayRead(localX,&xarr);CHKERRQ(ierr); 114*c4762a1bSJed Brown ierr = VecGetArrayRead(localXold,&xoldarr);CHKERRQ(ierr); 115*c4762a1bSJed Brown ierr = VecGetArray(localF,&farr);CHKERRQ(ierr); 116*c4762a1bSJed Brown 117*c4762a1bSJed Brown for (i=0; i<nv; i++) { 118*c4762a1bSJed Brown ierr = DMNetworkIsGhostVertex(networkdm,vtx[i],&ghostvtex);CHKERRQ(ierr); 119*c4762a1bSJed Brown if(ghostvtex) continue; 120*c4762a1bSJed Brown 121*c4762a1bSJed Brown ierr = DMNetworkGetVariableOffset(networkdm,vtx[i],&offset);CHKERRQ(ierr); 122*c4762a1bSJed Brown ierr = DMNetworkGetNumVariables(networkdm,vtx[i],&nvar);CHKERRQ(ierr); 123*c4762a1bSJed Brown for (j=0; j<nvar; j++) { 124*c4762a1bSJed Brown farr[offset+j] = xarr[offset+j] - xoldarr[offset+j]; 125*c4762a1bSJed Brown } 126*c4762a1bSJed Brown } 127*c4762a1bSJed Brown 128*c4762a1bSJed Brown ierr = VecRestoreArrayRead(localX,&xarr);CHKERRQ(ierr); 129*c4762a1bSJed Brown ierr = VecRestoreArrayRead(localXold,&xoldarr);CHKERRQ(ierr); 130*c4762a1bSJed Brown ierr = VecRestoreArray(localF,&farr);CHKERRQ(ierr); 131*c4762a1bSJed Brown PetscFunctionReturn(0); 132*c4762a1bSJed Brown } 133*c4762a1bSJed Brown 134*c4762a1bSJed Brown PetscErrorCode FormFunction(SNES snes,Vec X,Vec F,void *appctx) 135*c4762a1bSJed Brown { 136*c4762a1bSJed Brown PetscErrorCode ierr; 137*c4762a1bSJed Brown DM networkdm; 138*c4762a1bSJed Brown Vec localX,localF; 139*c4762a1bSJed Brown PetscInt nv,ne; 140*c4762a1bSJed Brown const PetscInt *vtx,*edges; 141*c4762a1bSJed Brown PetscMPIInt rank; 142*c4762a1bSJed Brown MPI_Comm comm; 143*c4762a1bSJed Brown UserCtx *user = (UserCtx*)appctx; 144*c4762a1bSJed Brown UserCtx_Power appctx_power = (*user).appctx_power; 145*c4762a1bSJed Brown 146*c4762a1bSJed Brown PetscFunctionBegin; 147*c4762a1bSJed Brown ierr = SNESGetDM(snes,&networkdm);CHKERRQ(ierr); 148*c4762a1bSJed Brown ierr = PetscObjectGetComm((PetscObject)networkdm,&comm);CHKERRQ(ierr); 149*c4762a1bSJed Brown ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 150*c4762a1bSJed Brown 151*c4762a1bSJed Brown ierr = DMGetLocalVector(networkdm,&localX);CHKERRQ(ierr); 152*c4762a1bSJed Brown ierr = DMGetLocalVector(networkdm,&localF);CHKERRQ(ierr); 153*c4762a1bSJed Brown ierr = VecSet(F,0.0);CHKERRQ(ierr); 154*c4762a1bSJed Brown ierr = VecSet(localF,0.0);CHKERRQ(ierr); 155*c4762a1bSJed Brown 156*c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(networkdm,X,INSERT_VALUES,localX);CHKERRQ(ierr); 157*c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(networkdm,X,INSERT_VALUES,localX);CHKERRQ(ierr); 158*c4762a1bSJed Brown 159*c4762a1bSJed Brown /* Form Function for power subnetwork */ 160*c4762a1bSJed Brown ierr = DMNetworkGetSubnetworkInfo(networkdm,0,&nv,&ne,&vtx,&edges);CHKERRQ(ierr); 161*c4762a1bSJed Brown if (user->subsnes_id == 1) { /* snes_water only */ 162*c4762a1bSJed Brown ierr = FormFunction_Dummy(networkdm,localX,localF,nv,ne,vtx,edges,user);CHKERRQ(ierr); 163*c4762a1bSJed Brown } else { 164*c4762a1bSJed Brown ierr = FormFunction_Power(networkdm,localX,localF,nv,ne,vtx,edges,&appctx_power);CHKERRQ(ierr); 165*c4762a1bSJed Brown } 166*c4762a1bSJed Brown 167*c4762a1bSJed Brown /* Form Function for water subnetwork */ 168*c4762a1bSJed Brown ierr = DMNetworkGetSubnetworkInfo(networkdm,1,&nv,&ne,&vtx,&edges);CHKERRQ(ierr); 169*c4762a1bSJed Brown if (user->subsnes_id == 0) { /* snes_power only */ 170*c4762a1bSJed Brown ierr = FormFunction_Dummy(networkdm,localX,localF,nv,ne,vtx,edges,user);CHKERRQ(ierr); 171*c4762a1bSJed Brown } else { 172*c4762a1bSJed Brown ierr = FormFunction_Water(networkdm,localX,localF,nv,ne,vtx,edges,NULL);CHKERRQ(ierr); 173*c4762a1bSJed Brown } 174*c4762a1bSJed Brown 175*c4762a1bSJed Brown #if 0 176*c4762a1bSJed Brown /* Form Function for the coupling subnetwork -- experimental, not done yet */ 177*c4762a1bSJed Brown ierr = DMNetworkGetSubnetworkCoupleInfo(networkdm,0,&ne,&edges);CHKERRQ(ierr); 178*c4762a1bSJed Brown if (ne) { 179*c4762a1bSJed Brown const PetscInt *cone,*connedges,*econe; 180*c4762a1bSJed Brown PetscInt key,vid[2],nconnedges,k,e,keye; 181*c4762a1bSJed Brown void* component; 182*c4762a1bSJed Brown AppCtx_Water appctx_water = (*user).appctx_water; 183*c4762a1bSJed Brown 184*c4762a1bSJed Brown ierr = DMNetworkGetConnectedVertices(networkdm,edges[0],&cone);CHKERRQ(ierr); 185*c4762a1bSJed Brown 186*c4762a1bSJed Brown ierr = DMNetworkGetGlobalVertexIndex(networkdm,cone[0],&vid[0]);CHKERRQ(ierr); 187*c4762a1bSJed Brown ierr = DMNetworkGetGlobalVertexIndex(networkdm,cone[1],&vid[1]);CHKERRQ(ierr); 188*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] Formfunction, coupling subnetwork: ne %d; connected vertices %d %d\n",rank,ne,vid[0],vid[1]);CHKERRQ(ierr); 189*c4762a1bSJed Brown 190*c4762a1bSJed Brown /* Get coupling powernet load vertex */ 191*c4762a1bSJed Brown ierr = DMNetworkGetComponent(networkdm,cone[0],1,&key,&component);CHKERRQ(ierr); 192*c4762a1bSJed Brown if (key != appctx_power.compkey_load) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Not a power load vertex"); 193*c4762a1bSJed Brown 194*c4762a1bSJed Brown /* Get coupling water vertex and pump edge */ 195*c4762a1bSJed Brown ierr = DMNetworkGetComponent(networkdm,cone[1],0,&key,&component);CHKERRQ(ierr); 196*c4762a1bSJed Brown if (key != appctx_water.compkey_vtx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Not a water vertex"); 197*c4762a1bSJed Brown 198*c4762a1bSJed Brown /* get its supporting edges */ 199*c4762a1bSJed Brown ierr = DMNetworkGetSupportingEdges(networkdm,cone[1],&nconnedges,&connedges);CHKERRQ(ierr); 200*c4762a1bSJed Brown 201*c4762a1bSJed Brown for (k=0; k<nconnedges; k++) { 202*c4762a1bSJed Brown e = connedges[k]; 203*c4762a1bSJed Brown ierr = DMNetworkGetComponent(networkdm,e,0,&keye,&component);CHKERRQ(ierr); 204*c4762a1bSJed Brown 205*c4762a1bSJed Brown if (keye == appctx_water.compkey_edge) { /* water_compkey_edge */ 206*c4762a1bSJed Brown EDGE_Water edge=(EDGE_Water)component; 207*c4762a1bSJed Brown if (edge->type == EDGE_TYPE_PUMP) { 208*c4762a1bSJed Brown /* compute flow_pump */ 209*c4762a1bSJed Brown PetscInt offsetnode1,offsetnode2,key_0,key_1; 210*c4762a1bSJed Brown const PetscScalar *xarr; 211*c4762a1bSJed Brown PetscScalar *farr; 212*c4762a1bSJed Brown VERTEX_Water vertexnode1,vertexnode2; 213*c4762a1bSJed Brown 214*c4762a1bSJed Brown ierr = DMNetworkGetConnectedVertices(networkdm,e,&econe);CHKERRQ(ierr); 215*c4762a1bSJed Brown ierr = DMNetworkGetGlobalVertexIndex(networkdm,econe[0],&vid[0]);CHKERRQ(ierr); 216*c4762a1bSJed Brown ierr = DMNetworkGetGlobalVertexIndex(networkdm,econe[1],&vid[1]);CHKERRQ(ierr); 217*c4762a1bSJed Brown 218*c4762a1bSJed Brown ierr = VecGetArray(localF,&farr);CHKERRQ(ierr); 219*c4762a1bSJed Brown ierr = DMNetworkGetVariableOffset(networkdm,econe[0],&offsetnode1);CHKERRQ(ierr); 220*c4762a1bSJed Brown ierr = DMNetworkGetVariableOffset(networkdm,econe[1],&offsetnode2);CHKERRQ(ierr); 221*c4762a1bSJed Brown 222*c4762a1bSJed Brown ierr = VecGetArrayRead(localX,&xarr);CHKERRQ(ierr); 223*c4762a1bSJed Brown #if 0 224*c4762a1bSJed Brown PetscScalar hf,ht; 225*c4762a1bSJed Brown Pump *pump; 226*c4762a1bSJed Brown pump = &edge->pump; 227*c4762a1bSJed Brown hf = xarr[offsetnode1]; 228*c4762a1bSJed Brown ht = xarr[offsetnode2]; 229*c4762a1bSJed Brown 230*c4762a1bSJed Brown PetscScalar flow = Flow_Pump(pump,hf,ht); 231*c4762a1bSJed Brown PetscScalar Hp = 0.1; /* load->pl */ 232*c4762a1bSJed Brown PetscScalar flow_couple = 8.81*Hp*1.e6/(ht-hf); /* pump->h0; */ 233*c4762a1bSJed Brown /* ierr = PetscPrintf(PETSC_COMM_SELF,"pump %d: connected vtx %d %d; flow_pump %g flow_couple %g; offset %d %d\n",e,vid[0],vid[1],flow,flow_couple,offsetnode1,offsetnode2);CHKERRQ(ierr); */ 234*c4762a1bSJed Brown #endif 235*c4762a1bSJed Brown /* Get the components at the two vertices */ 236*c4762a1bSJed Brown ierr = DMNetworkGetComponent(networkdm,econe[0],0,&key_0,(void**)&vertexnode1);CHKERRQ(ierr); 237*c4762a1bSJed Brown if (key_0 != appctx_water.compkey_vtx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Not a water vertex"); 238*c4762a1bSJed Brown ierr = DMNetworkGetComponent(networkdm,econe[1],0,&key_1,(void**)&vertexnode2);CHKERRQ(ierr); 239*c4762a1bSJed Brown if (key_1 != appctx_water.compkey_vtx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Not a water vertex"); 240*c4762a1bSJed Brown #if 0 241*c4762a1bSJed Brown /* subtract flow_pump computed in FormFunction_Water() and add flow_couple to connected nodes */ 242*c4762a1bSJed Brown if (vertexnode1->type == VERTEX_TYPE_JUNCTION) { 243*c4762a1bSJed Brown farr[offsetnode1] += flow; 244*c4762a1bSJed Brown farr[offsetnode1] -= flow_couple; 245*c4762a1bSJed Brown } 246*c4762a1bSJed Brown if (vertexnode2->type == VERTEX_TYPE_JUNCTION) { 247*c4762a1bSJed Brown farr[offsetnode2] -= flow; 248*c4762a1bSJed Brown farr[offsetnode2] += flow_couple; 249*c4762a1bSJed Brown } 250*c4762a1bSJed Brown #endif 251*c4762a1bSJed Brown ierr = VecRestoreArrayRead(localX,&xarr);CHKERRQ(ierr); 252*c4762a1bSJed Brown ierr = VecRestoreArray(localF,&farr);CHKERRQ(ierr); 253*c4762a1bSJed Brown } 254*c4762a1bSJed Brown } 255*c4762a1bSJed Brown } 256*c4762a1bSJed Brown } 257*c4762a1bSJed Brown #endif 258*c4762a1bSJed Brown 259*c4762a1bSJed Brown ierr = DMRestoreLocalVector(networkdm,&localX);CHKERRQ(ierr); 260*c4762a1bSJed Brown 261*c4762a1bSJed Brown ierr = DMLocalToGlobalBegin(networkdm,localF,ADD_VALUES,F);CHKERRQ(ierr); 262*c4762a1bSJed Brown ierr = DMLocalToGlobalEnd(networkdm,localF,ADD_VALUES,F);CHKERRQ(ierr); 263*c4762a1bSJed Brown ierr = DMRestoreLocalVector(networkdm,&localF);CHKERRQ(ierr); 264*c4762a1bSJed Brown /* ierr = VecView(F,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */ 265*c4762a1bSJed Brown PetscFunctionReturn(0); 266*c4762a1bSJed Brown } 267*c4762a1bSJed Brown 268*c4762a1bSJed Brown PetscErrorCode SetInitialGuess(DM networkdm,Vec X,void* appctx) 269*c4762a1bSJed Brown { 270*c4762a1bSJed Brown PetscErrorCode ierr; 271*c4762a1bSJed Brown PetscInt nv,ne; 272*c4762a1bSJed Brown const PetscInt *vtx,*edges; 273*c4762a1bSJed Brown UserCtx *user = (UserCtx*)appctx; 274*c4762a1bSJed Brown Vec localX = user->localXold; 275*c4762a1bSJed Brown UserCtx_Power appctx_power = (*user).appctx_power; 276*c4762a1bSJed Brown 277*c4762a1bSJed Brown PetscFunctionBegin; 278*c4762a1bSJed Brown ierr = VecSet(X,0.0);CHKERRQ(ierr); 279*c4762a1bSJed Brown ierr = VecSet(localX,0.0);CHKERRQ(ierr); 280*c4762a1bSJed Brown 281*c4762a1bSJed Brown /* Set initial guess for power subnetwork */ 282*c4762a1bSJed Brown ierr = DMNetworkGetSubnetworkInfo(networkdm,0,&nv,&ne,&vtx,&edges);CHKERRQ(ierr); 283*c4762a1bSJed Brown ierr = SetInitialGuess_Power(networkdm,localX,nv,ne,vtx,edges,&appctx_power);CHKERRQ(ierr); 284*c4762a1bSJed Brown 285*c4762a1bSJed Brown /* Set initial guess for water subnetwork */ 286*c4762a1bSJed Brown ierr = DMNetworkGetSubnetworkInfo(networkdm,1,&nv,&ne,&vtx,&edges);CHKERRQ(ierr); 287*c4762a1bSJed Brown ierr = SetInitialGuess_Water(networkdm,localX,nv,ne,vtx,edges,NULL);CHKERRQ(ierr); 288*c4762a1bSJed Brown 289*c4762a1bSJed Brown #if 0 290*c4762a1bSJed Brown /* Set initial guess for the coupling subnet */ 291*c4762a1bSJed Brown ierr = DMNetworkGetSubnetworkInfo(networkdm,2,&nv,&ne,&vtx,&edges);CHKERRQ(ierr); 292*c4762a1bSJed Brown if (ne) { 293*c4762a1bSJed Brown const PetscInt *cone; 294*c4762a1bSJed Brown ierr = DMNetworkGetConnectedVertices(networkdm,edges[0],&cone);CHKERRQ(ierr); 295*c4762a1bSJed Brown } 296*c4762a1bSJed Brown #endif 297*c4762a1bSJed Brown 298*c4762a1bSJed Brown ierr = DMLocalToGlobalBegin(networkdm,localX,ADD_VALUES,X);CHKERRQ(ierr); 299*c4762a1bSJed Brown ierr = DMLocalToGlobalEnd(networkdm,localX,ADD_VALUES,X);CHKERRQ(ierr); 300*c4762a1bSJed Brown PetscFunctionReturn(0); 301*c4762a1bSJed Brown } 302*c4762a1bSJed Brown 303*c4762a1bSJed Brown int main(int argc,char **argv) 304*c4762a1bSJed Brown { 305*c4762a1bSJed Brown PetscErrorCode ierr; 306*c4762a1bSJed Brown DM networkdm; 307*c4762a1bSJed Brown PetscLogStage stage[4]; 308*c4762a1bSJed Brown PetscMPIInt rank; 309*c4762a1bSJed Brown PetscInt nsubnet=2,nsubnetCouple=1,numVertices[2],numEdges[2],numEdgesCouple[1]; 310*c4762a1bSJed Brown PetscInt i,j,nv,ne; 311*c4762a1bSJed Brown PetscInt *edgelist[2]; 312*c4762a1bSJed Brown const PetscInt *vtx,*edges; 313*c4762a1bSJed Brown Vec X,F; 314*c4762a1bSJed Brown SNES snes,snes_power,snes_water; 315*c4762a1bSJed Brown Mat Jac; 316*c4762a1bSJed Brown PetscBool viewJ=PETSC_FALSE,viewX=PETSC_FALSE,viewDM=PETSC_FALSE,test=PETSC_FALSE; 317*c4762a1bSJed Brown UserCtx user; 318*c4762a1bSJed Brown PetscInt it_max=10; 319*c4762a1bSJed Brown SNESConvergedReason reason; 320*c4762a1bSJed Brown 321*c4762a1bSJed Brown /* Power subnetwork */ 322*c4762a1bSJed Brown UserCtx_Power *appctx_power = &user.appctx_power; 323*c4762a1bSJed Brown char pfdata_file[PETSC_MAX_PATH_LEN]="power/case9.m"; 324*c4762a1bSJed Brown PFDATA *pfdata=NULL; 325*c4762a1bSJed Brown PetscInt genj,loadj; 326*c4762a1bSJed Brown PetscInt *edgelist_power=NULL; 327*c4762a1bSJed Brown PetscScalar Sbase=0.0; 328*c4762a1bSJed Brown 329*c4762a1bSJed Brown /* Water subnetwork */ 330*c4762a1bSJed Brown AppCtx_Water *appctx_water = &user.appctx_water; 331*c4762a1bSJed Brown WATERDATA *waterdata=NULL; 332*c4762a1bSJed Brown char waterdata_file[PETSC_MAX_PATH_LEN]="water/sample1.inp"; 333*c4762a1bSJed Brown PetscInt *edgelist_water=NULL; 334*c4762a1bSJed Brown 335*c4762a1bSJed Brown /* Coupling subnetwork */ 336*c4762a1bSJed Brown PetscInt *edgelist_couple=NULL; 337*c4762a1bSJed Brown 338*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,"ex1options",help);if (ierr) return ierr; 339*c4762a1bSJed Brown ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr); 340*c4762a1bSJed Brown 341*c4762a1bSJed Brown /* (1) Read Data - Only rank 0 reads the data */ 342*c4762a1bSJed Brown /*--------------------------------------------*/ 343*c4762a1bSJed Brown ierr = PetscLogStageRegister("Read Data",&stage[0]);CHKERRQ(ierr); 344*c4762a1bSJed Brown PetscLogStagePush(stage[0]); 345*c4762a1bSJed Brown 346*c4762a1bSJed Brown for (i=0; i<nsubnet; i++) { 347*c4762a1bSJed Brown numVertices[i] = 0; 348*c4762a1bSJed Brown numEdges[i] = 0; 349*c4762a1bSJed Brown } 350*c4762a1bSJed Brown for (i=0; i<nsubnetCouple; i++) { 351*c4762a1bSJed Brown numEdgesCouple[0] = 0; 352*c4762a1bSJed Brown } 353*c4762a1bSJed Brown 354*c4762a1bSJed Brown /* READ THE DATA FOR THE FIRST SUBNETWORK: Electric Power Grid */ 355*c4762a1bSJed Brown if (!rank) { 356*c4762a1bSJed Brown ierr = PetscOptionsGetString(NULL,NULL,"-pfdata",pfdata_file,PETSC_MAX_PATH_LEN-1,NULL);CHKERRQ(ierr); 357*c4762a1bSJed Brown ierr = PetscNew(&pfdata);CHKERRQ(ierr); 358*c4762a1bSJed Brown ierr = PFReadMatPowerData(pfdata,pfdata_file);CHKERRQ(ierr); 359*c4762a1bSJed Brown Sbase = pfdata->sbase; 360*c4762a1bSJed Brown 361*c4762a1bSJed Brown numEdges[0] = pfdata->nbranch; 362*c4762a1bSJed Brown numVertices[0] = pfdata->nbus; 363*c4762a1bSJed Brown 364*c4762a1bSJed Brown ierr = PetscMalloc1(2*numEdges[0],&edgelist_power);CHKERRQ(ierr); 365*c4762a1bSJed Brown ierr = GetListofEdges_Power(pfdata,edgelist_power);CHKERRQ(ierr); 366*c4762a1bSJed Brown #if 0 367*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"edgelist_power:\n");CHKERRQ(ierr); 368*c4762a1bSJed Brown for (i=0; i<numEdges[0]; i++) { 369*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"[%D %D]",edgelist_power[2*i],edgelist_power[2*i+1]);CHKERRQ(ierr); 370*c4762a1bSJed Brown } 371*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"\n");CHKERRQ(ierr); 372*c4762a1bSJed Brown #endif 373*c4762a1bSJed Brown } 374*c4762a1bSJed Brown /* Broadcast power Sbase to all processors */ 375*c4762a1bSJed Brown ierr = MPI_Bcast(&Sbase,1,MPIU_SCALAR,0,PETSC_COMM_WORLD);CHKERRQ(ierr); 376*c4762a1bSJed Brown appctx_power->Sbase = Sbase; 377*c4762a1bSJed Brown appctx_power->jac_error = PETSC_FALSE; 378*c4762a1bSJed Brown /* If external option activated. Introduce error in jacobian */ 379*c4762a1bSJed Brown ierr = PetscOptionsHasName(NULL,NULL, "-jac_error", &appctx_power->jac_error);CHKERRQ(ierr); 380*c4762a1bSJed Brown 381*c4762a1bSJed Brown /* GET DATA FOR THE SECOND SUBNETWORK: Water */ 382*c4762a1bSJed Brown if (!rank) { 383*c4762a1bSJed Brown ierr = PetscNew(&waterdata);CHKERRQ(ierr); 384*c4762a1bSJed Brown ierr = PetscOptionsGetString(NULL,NULL,"-waterdata",waterdata_file,PETSC_MAX_PATH_LEN-1,NULL);CHKERRQ(ierr); 385*c4762a1bSJed Brown ierr = WaterReadData(waterdata,waterdata_file);CHKERRQ(ierr); 386*c4762a1bSJed Brown 387*c4762a1bSJed Brown ierr = PetscCalloc1(2*waterdata->nedge,&edgelist_water);CHKERRQ(ierr); 388*c4762a1bSJed Brown ierr = GetListofEdges_Water(waterdata,edgelist_water);CHKERRQ(ierr); 389*c4762a1bSJed Brown numEdges[1] = waterdata->nedge; 390*c4762a1bSJed Brown numVertices[1] = waterdata->nvertex; 391*c4762a1bSJed Brown #if 0 392*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"edgelist_water:\n");CHKERRQ(ierr); 393*c4762a1bSJed Brown for (i=0; i<numEdges[1]; i++) { 394*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"[%D %D]",edgelist_water[2*i],edgelist_water[2*i+1]);CHKERRQ(ierr); 395*c4762a1bSJed Brown } 396*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,("\n");CHKERRQ(ierr); 397*c4762a1bSJed Brown #endif 398*c4762a1bSJed Brown } 399*c4762a1bSJed Brown 400*c4762a1bSJed Brown /* Get data for the coupling subnetwork */ 401*c4762a1bSJed Brown if (!rank) { 402*c4762a1bSJed Brown numEdgesCouple[0] = 1; 403*c4762a1bSJed Brown 404*c4762a1bSJed Brown ierr = PetscMalloc1(4*numEdgesCouple[0],&edgelist_couple);CHKERRQ(ierr); 405*c4762a1bSJed Brown edgelist_couple[0] = 0; edgelist_couple[1] = 4; /* from node: net[0] vertex[4] */ 406*c4762a1bSJed Brown edgelist_couple[2] = 1; edgelist_couple[3] = 0; /* to node: net[1] vertex[0] */ 407*c4762a1bSJed Brown } 408*c4762a1bSJed Brown PetscLogStagePop(); 409*c4762a1bSJed Brown 410*c4762a1bSJed Brown /* (2) Create network */ 411*c4762a1bSJed Brown /*--------------------*/ 412*c4762a1bSJed Brown ierr = MPI_Barrier(PETSC_COMM_WORLD);CHKERRQ(ierr); 413*c4762a1bSJed Brown ierr = PetscLogStageRegister("Net Setup",&stage[1]);CHKERRQ(ierr); 414*c4762a1bSJed Brown PetscLogStagePush(stage[1]); 415*c4762a1bSJed Brown 416*c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-viewDM",&viewDM,NULL);CHKERRQ(ierr); 417*c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-test",&test,NULL);CHKERRQ(ierr); 418*c4762a1bSJed Brown 419*c4762a1bSJed Brown /* Create an empty network object */ 420*c4762a1bSJed Brown ierr = DMNetworkCreate(PETSC_COMM_WORLD,&networkdm);CHKERRQ(ierr); 421*c4762a1bSJed Brown 422*c4762a1bSJed Brown /* Register the components in the network */ 423*c4762a1bSJed Brown ierr = DMNetworkRegisterComponent(networkdm,"branchstruct",sizeof(struct _p_EDGE_Power),&appctx_power->compkey_branch);CHKERRQ(ierr); 424*c4762a1bSJed Brown ierr = DMNetworkRegisterComponent(networkdm,"busstruct",sizeof(struct _p_VERTEX_Power),&appctx_power->compkey_bus);CHKERRQ(ierr); 425*c4762a1bSJed Brown ierr = DMNetworkRegisterComponent(networkdm,"genstruct",sizeof(struct _p_GEN),&appctx_power->compkey_gen);CHKERRQ(ierr); 426*c4762a1bSJed Brown ierr = DMNetworkRegisterComponent(networkdm,"loadstruct",sizeof(struct _p_LOAD),&appctx_power->compkey_load);CHKERRQ(ierr); 427*c4762a1bSJed Brown 428*c4762a1bSJed Brown ierr = DMNetworkRegisterComponent(networkdm,"edge_water",sizeof(struct _p_EDGE_Water),&appctx_water->compkey_edge);CHKERRQ(ierr); 429*c4762a1bSJed Brown ierr = DMNetworkRegisterComponent(networkdm,"vertex_water",sizeof(struct _p_VERTEX_Water),&appctx_water->compkey_vtx);CHKERRQ(ierr); 430*c4762a1bSJed Brown 431*c4762a1bSJed Brown if (!rank) { 432*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] Total local nvertices %D + %D = %D, nedges %D + %D + %D = %D\n",rank,numVertices[0],numVertices[1],numVertices[0]+numVertices[1],numEdges[0],numEdges[1],numEdgesCouple[0],numEdges[0]+numEdges[1]+numEdgesCouple[0]);CHKERRQ(ierr); 433*c4762a1bSJed Brown } 434*c4762a1bSJed Brown 435*c4762a1bSJed Brown ierr = DMNetworkSetSizes(networkdm,nsubnet,numVertices,numEdges,nsubnetCouple,numEdgesCouple);CHKERRQ(ierr); 436*c4762a1bSJed Brown 437*c4762a1bSJed Brown /* Add edge connectivity */ 438*c4762a1bSJed Brown edgelist[0] = edgelist_power; 439*c4762a1bSJed Brown edgelist[1] = edgelist_water; 440*c4762a1bSJed Brown ierr = DMNetworkSetEdgeList(networkdm,edgelist,&edgelist_couple);CHKERRQ(ierr); 441*c4762a1bSJed Brown 442*c4762a1bSJed Brown /* Set up the network layout */ 443*c4762a1bSJed Brown ierr = DMNetworkLayoutSetUp(networkdm);CHKERRQ(ierr); 444*c4762a1bSJed Brown 445*c4762a1bSJed Brown /* Add network components - only process[0] has any data to add */ 446*c4762a1bSJed Brown /* ADD VARIABLES AND COMPONENTS FOR THE POWER SUBNETWORK */ 447*c4762a1bSJed Brown genj = 0; loadj = 0; 448*c4762a1bSJed Brown if (!rank) { 449*c4762a1bSJed Brown ierr = DMNetworkGetSubnetworkInfo(networkdm,0,&nv,&ne,&vtx,&edges);CHKERRQ(ierr); 450*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Power network: nv %D, ne %D\n",nv,ne);CHKERRQ(ierr); 451*c4762a1bSJed Brown for (i = 0; i < ne; i++) { 452*c4762a1bSJed Brown ierr = DMNetworkAddComponent(networkdm,edges[i],appctx_power->compkey_branch,&pfdata->branch[i]);CHKERRQ(ierr); 453*c4762a1bSJed Brown } 454*c4762a1bSJed Brown 455*c4762a1bSJed Brown for (i = 0; i < nv; i++) { 456*c4762a1bSJed Brown ierr = DMNetworkAddComponent(networkdm,vtx[i],appctx_power->compkey_bus,&pfdata->bus[i]);CHKERRQ(ierr); 457*c4762a1bSJed Brown if (pfdata->bus[i].ngen) { 458*c4762a1bSJed Brown for (j = 0; j < pfdata->bus[i].ngen; j++) { 459*c4762a1bSJed Brown ierr = DMNetworkAddComponent(networkdm,vtx[i],appctx_power->compkey_gen,&pfdata->gen[genj++]);CHKERRQ(ierr); 460*c4762a1bSJed Brown } 461*c4762a1bSJed Brown } 462*c4762a1bSJed Brown if (pfdata->bus[i].nload) { 463*c4762a1bSJed Brown for (j=0; j < pfdata->bus[i].nload; j++) { 464*c4762a1bSJed Brown ierr = DMNetworkAddComponent(networkdm,vtx[i],appctx_power->compkey_load,&pfdata->load[loadj++]);CHKERRQ(ierr); 465*c4762a1bSJed Brown } 466*c4762a1bSJed Brown } 467*c4762a1bSJed Brown /* Add number of variables */ 468*c4762a1bSJed Brown ierr = DMNetworkAddNumVariables(networkdm,vtx[i],2);CHKERRQ(ierr); 469*c4762a1bSJed Brown } 470*c4762a1bSJed Brown } 471*c4762a1bSJed Brown 472*c4762a1bSJed Brown /* ADD VARIABLES AND COMPONENTS FOR THE WATER SUBNETWORK */ 473*c4762a1bSJed Brown if (!rank) { 474*c4762a1bSJed Brown ierr = DMNetworkGetSubnetworkInfo(networkdm,1,&nv,&ne,&vtx,&edges);CHKERRQ(ierr); 475*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Water network: nv %D, ne %D\n",nv,ne);CHKERRQ(ierr); 476*c4762a1bSJed Brown for (i = 0; i < ne; i++) { 477*c4762a1bSJed Brown ierr = DMNetworkAddComponent(networkdm,edges[i],appctx_water->compkey_edge,&waterdata->edge[i]);CHKERRQ(ierr); 478*c4762a1bSJed Brown } 479*c4762a1bSJed Brown 480*c4762a1bSJed Brown for (i = 0; i < nv; i++) { 481*c4762a1bSJed Brown ierr = DMNetworkAddComponent(networkdm,vtx[i],appctx_water->compkey_vtx,&waterdata->vertex[i]);CHKERRQ(ierr); 482*c4762a1bSJed Brown /* Add number of variables */ 483*c4762a1bSJed Brown ierr = DMNetworkAddNumVariables(networkdm,vtx[i],1);CHKERRQ(ierr); 484*c4762a1bSJed Brown } 485*c4762a1bSJed Brown } 486*c4762a1bSJed Brown 487*c4762a1bSJed Brown /* Set up DM for use */ 488*c4762a1bSJed Brown ierr = DMSetUp(networkdm);CHKERRQ(ierr); 489*c4762a1bSJed Brown if (viewDM) { 490*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"\nAfter DMSetUp, DMView:\n");CHKERRQ(ierr); 491*c4762a1bSJed Brown ierr = DMView(networkdm,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 492*c4762a1bSJed Brown } 493*c4762a1bSJed Brown 494*c4762a1bSJed Brown ierr = DMNetworkGetSubnetworkCoupleInfo(networkdm,0,&ne,&edges);CHKERRQ(ierr); 495*c4762a1bSJed Brown if (ne && viewDM) { 496*c4762a1bSJed Brown const PetscInt *cone; 497*c4762a1bSJed Brown PetscInt vid[2]; 498*c4762a1bSJed Brown ierr = DMNetworkGetConnectedVertices(networkdm,edges[0],&cone);CHKERRQ(ierr); 499*c4762a1bSJed Brown 500*c4762a1bSJed Brown ierr = DMNetworkGetGlobalVertexIndex(networkdm,cone[0],&vid[0]);CHKERRQ(ierr); 501*c4762a1bSJed Brown ierr = DMNetworkGetGlobalVertexIndex(networkdm,cone[1],&vid[1]);CHKERRQ(ierr); 502*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] Coupling subnetwork: ne %d; connected vertices %d %d\n",rank,ne,vid[0],vid[1]);CHKERRQ(ierr); 503*c4762a1bSJed Brown } 504*c4762a1bSJed Brown 505*c4762a1bSJed Brown /* Free user objects */ 506*c4762a1bSJed Brown if (!rank) { 507*c4762a1bSJed Brown ierr = PetscFree(edgelist_power);CHKERRQ(ierr); 508*c4762a1bSJed Brown ierr = PetscFree(pfdata->bus);CHKERRQ(ierr); 509*c4762a1bSJed Brown ierr = PetscFree(pfdata->gen);CHKERRQ(ierr); 510*c4762a1bSJed Brown ierr = PetscFree(pfdata->branch);CHKERRQ(ierr); 511*c4762a1bSJed Brown ierr = PetscFree(pfdata->load);CHKERRQ(ierr); 512*c4762a1bSJed Brown ierr = PetscFree(pfdata);CHKERRQ(ierr); 513*c4762a1bSJed Brown 514*c4762a1bSJed Brown ierr = PetscFree(edgelist_water);CHKERRQ(ierr); 515*c4762a1bSJed Brown ierr = PetscFree(waterdata->vertex);CHKERRQ(ierr); 516*c4762a1bSJed Brown ierr = PetscFree(waterdata->edge);CHKERRQ(ierr); 517*c4762a1bSJed Brown 518*c4762a1bSJed Brown ierr = PetscFree(edgelist_couple);CHKERRQ(ierr); 519*c4762a1bSJed Brown ierr = PetscFree(waterdata);CHKERRQ(ierr); 520*c4762a1bSJed Brown } 521*c4762a1bSJed Brown 522*c4762a1bSJed Brown /* Distribute networkdm to multiple processes */ 523*c4762a1bSJed Brown ierr = DMNetworkDistribute(&networkdm,0);CHKERRQ(ierr); 524*c4762a1bSJed Brown if (viewDM) { 525*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"\nAfter DMNetworkDistribute, DMView:\n");CHKERRQ(ierr); 526*c4762a1bSJed Brown ierr = DMView(networkdm,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 527*c4762a1bSJed Brown } 528*c4762a1bSJed Brown 529*c4762a1bSJed Brown /* Test DMNetworkGetSubnetworkCoupleInfo() */ 530*c4762a1bSJed Brown ierr = MPI_Barrier(PETSC_COMM_WORLD);CHKERRQ(ierr); 531*c4762a1bSJed Brown if (test) { 532*c4762a1bSJed Brown for (i=0; i<nsubnetCouple; i++) { 533*c4762a1bSJed Brown ierr = DMNetworkGetSubnetworkCoupleInfo(networkdm,0,&ne,&edges);CHKERRQ(ierr); 534*c4762a1bSJed Brown if (ne && viewDM) { 535*c4762a1bSJed Brown const PetscInt *cone; 536*c4762a1bSJed Brown PetscInt vid[2]; 537*c4762a1bSJed Brown ierr = DMNetworkGetConnectedVertices(networkdm,edges[0],&cone);CHKERRQ(ierr); 538*c4762a1bSJed Brown 539*c4762a1bSJed Brown ierr = DMNetworkGetGlobalVertexIndex(networkdm,cone[0],&vid[0]);CHKERRQ(ierr); 540*c4762a1bSJed Brown ierr = DMNetworkGetGlobalVertexIndex(networkdm,cone[1],&vid[1]);CHKERRQ(ierr); 541*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] After DMNetworkDistribute(), subnetworkCouple[%d]: ne %d; connected vertices %d %d\n",rank,i,ne,vid[0],vid[1]);CHKERRQ(ierr); 542*c4762a1bSJed Brown } 543*c4762a1bSJed Brown } 544*c4762a1bSJed Brown } 545*c4762a1bSJed Brown 546*c4762a1bSJed Brown ierr = DMCreateGlobalVector(networkdm,&X);CHKERRQ(ierr); 547*c4762a1bSJed Brown ierr = VecDuplicate(X,&F);CHKERRQ(ierr); 548*c4762a1bSJed Brown ierr = DMGetLocalVector(networkdm,&user.localXold);CHKERRQ(ierr); 549*c4762a1bSJed Brown 550*c4762a1bSJed Brown PetscLogStagePop(); 551*c4762a1bSJed Brown 552*c4762a1bSJed Brown /* (3) Setup Solvers */ 553*c4762a1bSJed Brown /*-------------------*/ 554*c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-viewJ",&viewJ,NULL);CHKERRQ(ierr); 555*c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-viewX",&viewX,NULL);CHKERRQ(ierr); 556*c4762a1bSJed Brown 557*c4762a1bSJed Brown ierr = PetscLogStageRegister("SNES Setup",&stage[2]);CHKERRQ(ierr); 558*c4762a1bSJed Brown PetscLogStagePush(stage[2]); 559*c4762a1bSJed Brown 560*c4762a1bSJed Brown ierr = SetInitialGuess(networkdm,X,&user);CHKERRQ(ierr); 561*c4762a1bSJed Brown /* ierr = VecView(X,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */ 562*c4762a1bSJed Brown 563*c4762a1bSJed Brown /* Create coupled snes */ 564*c4762a1bSJed Brown /*-------------------- */ 565*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"SNES_coupled setup ......\n");CHKERRQ(ierr); 566*c4762a1bSJed Brown user.subsnes_id = nsubnet; 567*c4762a1bSJed Brown ierr = SNESCreate(PETSC_COMM_WORLD,&snes);CHKERRQ(ierr); 568*c4762a1bSJed Brown ierr = SNESSetDM(snes,networkdm);CHKERRQ(ierr); 569*c4762a1bSJed Brown ierr = SNESSetOptionsPrefix(snes,"coupled_");CHKERRQ(ierr); 570*c4762a1bSJed Brown ierr = SNESSetFunction(snes,F,FormFunction,&user);CHKERRQ(ierr); 571*c4762a1bSJed Brown ierr = SNESMonitorSet(snes,UserMonitor,&user,NULL);CHKERRQ(ierr); 572*c4762a1bSJed Brown ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 573*c4762a1bSJed Brown 574*c4762a1bSJed Brown if (viewJ) { 575*c4762a1bSJed Brown /* View Jac structure */ 576*c4762a1bSJed Brown ierr = SNESGetJacobian(snes,&Jac,NULL,NULL,NULL);CHKERRQ(ierr); 577*c4762a1bSJed Brown ierr = MatView(Jac,PETSC_VIEWER_DRAW_WORLD);CHKERRQ(ierr); 578*c4762a1bSJed Brown } 579*c4762a1bSJed Brown 580*c4762a1bSJed Brown if (viewX) { 581*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Solution:\n");CHKERRQ(ierr); 582*c4762a1bSJed Brown ierr = VecView(X,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 583*c4762a1bSJed Brown } 584*c4762a1bSJed Brown 585*c4762a1bSJed Brown if (viewJ) { 586*c4762a1bSJed Brown /* View assembled Jac */ 587*c4762a1bSJed Brown ierr = MatView(Jac,PETSC_VIEWER_DRAW_WORLD);CHKERRQ(ierr); 588*c4762a1bSJed Brown } 589*c4762a1bSJed Brown 590*c4762a1bSJed Brown /* Create snes_power */ 591*c4762a1bSJed Brown /*-------------------*/ 592*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"SNES_power setup ......\n");CHKERRQ(ierr); 593*c4762a1bSJed Brown ierr = SetInitialGuess(networkdm,X,&user);CHKERRQ(ierr); 594*c4762a1bSJed Brown 595*c4762a1bSJed Brown user.subsnes_id = 0; 596*c4762a1bSJed Brown ierr = SNESCreate(PETSC_COMM_WORLD,&snes_power);CHKERRQ(ierr); 597*c4762a1bSJed Brown ierr = SNESSetDM(snes_power,networkdm);CHKERRQ(ierr); 598*c4762a1bSJed Brown ierr = SNESSetOptionsPrefix(snes_power,"power_");CHKERRQ(ierr); 599*c4762a1bSJed Brown ierr = SNESSetFunction(snes_power,F,FormFunction,&user);CHKERRQ(ierr); 600*c4762a1bSJed Brown ierr = SNESMonitorSet(snes_power,UserMonitor,&user,NULL);CHKERRQ(ierr); 601*c4762a1bSJed Brown 602*c4762a1bSJed Brown /* Use user-provide Jacobian */ 603*c4762a1bSJed Brown ierr = DMCreateMatrix(networkdm,&Jac);CHKERRQ(ierr); 604*c4762a1bSJed Brown ierr = SNESSetJacobian(snes_power,Jac,Jac,FormJacobian_subPower,&user);CHKERRQ(ierr); 605*c4762a1bSJed Brown ierr = SNESSetFromOptions(snes_power);CHKERRQ(ierr); 606*c4762a1bSJed Brown 607*c4762a1bSJed Brown if (viewX) { 608*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Power Solution:\n");CHKERRQ(ierr); 609*c4762a1bSJed Brown ierr = VecView(X,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 610*c4762a1bSJed Brown } 611*c4762a1bSJed Brown if (viewJ) { 612*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Power Jac:\n");CHKERRQ(ierr); 613*c4762a1bSJed Brown ierr = SNESGetJacobian(snes_power,&Jac,NULL,NULL,NULL);CHKERRQ(ierr); 614*c4762a1bSJed Brown ierr = MatView(Jac,PETSC_VIEWER_DRAW_WORLD);CHKERRQ(ierr); 615*c4762a1bSJed Brown /* ierr = MatView(Jac,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */ 616*c4762a1bSJed Brown } 617*c4762a1bSJed Brown 618*c4762a1bSJed Brown /* Create snes_water */ 619*c4762a1bSJed Brown /*-------------------*/ 620*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"SNES_water setup......\n");CHKERRQ(ierr); 621*c4762a1bSJed Brown ierr = SetInitialGuess(networkdm,X,&user);CHKERRQ(ierr); 622*c4762a1bSJed Brown 623*c4762a1bSJed Brown user.subsnes_id = 1; 624*c4762a1bSJed Brown ierr = SNESCreate(PETSC_COMM_WORLD,&snes_water);CHKERRQ(ierr); 625*c4762a1bSJed Brown ierr = SNESSetDM(snes_water,networkdm);CHKERRQ(ierr); 626*c4762a1bSJed Brown ierr = SNESSetOptionsPrefix(snes_water,"water_");CHKERRQ(ierr); 627*c4762a1bSJed Brown ierr = SNESSetFunction(snes_water,F,FormFunction,&user);CHKERRQ(ierr); 628*c4762a1bSJed Brown ierr = SNESMonitorSet(snes_water,UserMonitor,&user,NULL);CHKERRQ(ierr); 629*c4762a1bSJed Brown ierr = SNESSetFromOptions(snes_water);CHKERRQ(ierr); 630*c4762a1bSJed Brown 631*c4762a1bSJed Brown if (viewX) { 632*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Water Solution:\n");CHKERRQ(ierr); 633*c4762a1bSJed Brown ierr = VecView(X,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 634*c4762a1bSJed Brown } 635*c4762a1bSJed Brown PetscLogStagePop(); 636*c4762a1bSJed Brown 637*c4762a1bSJed Brown /* (4) Solve */ 638*c4762a1bSJed Brown /*-----------*/ 639*c4762a1bSJed Brown ierr = PetscLogStageRegister("SNES Solve",&stage[3]);CHKERRQ(ierr); 640*c4762a1bSJed Brown PetscLogStagePush(stage[3]); 641*c4762a1bSJed Brown user.it = 0; 642*c4762a1bSJed Brown reason = SNES_DIVERGED_DTOL; 643*c4762a1bSJed Brown while (user.it < it_max && (PetscInt)reason<0) { 644*c4762a1bSJed Brown #if 0 645*c4762a1bSJed Brown user.subsnes_id = 0; 646*c4762a1bSJed Brown ierr = SNESSolve(snes_power,NULL,X);CHKERRQ(ierr); 647*c4762a1bSJed Brown 648*c4762a1bSJed Brown user.subsnes_id = 1; 649*c4762a1bSJed Brown ierr = SNESSolve(snes_water,NULL,X);CHKERRQ(ierr); 650*c4762a1bSJed Brown #endif 651*c4762a1bSJed Brown user.subsnes_id = nsubnet; 652*c4762a1bSJed Brown ierr = SNESSolve(snes,NULL,X);CHKERRQ(ierr); 653*c4762a1bSJed Brown 654*c4762a1bSJed Brown ierr = SNESGetConvergedReason(snes,&reason);CHKERRQ(ierr); 655*c4762a1bSJed Brown user.it++; 656*c4762a1bSJed Brown } 657*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Coupled_SNES converged in %D iterations\n",user.it);CHKERRQ(ierr); 658*c4762a1bSJed Brown if (viewX) { 659*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Final Solution:\n");CHKERRQ(ierr); 660*c4762a1bSJed Brown ierr = VecView(X,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 661*c4762a1bSJed Brown } 662*c4762a1bSJed Brown PetscLogStagePop(); 663*c4762a1bSJed Brown 664*c4762a1bSJed Brown /* Free objects */ 665*c4762a1bSJed Brown /* -------------*/ 666*c4762a1bSJed Brown ierr = VecDestroy(&X);CHKERRQ(ierr); 667*c4762a1bSJed Brown ierr = VecDestroy(&F);CHKERRQ(ierr); 668*c4762a1bSJed Brown ierr = DMRestoreLocalVector(networkdm,&user.localXold);CHKERRQ(ierr); 669*c4762a1bSJed Brown 670*c4762a1bSJed Brown ierr = SNESDestroy(&snes);CHKERRQ(ierr); 671*c4762a1bSJed Brown ierr = MatDestroy(&Jac);CHKERRQ(ierr); 672*c4762a1bSJed Brown ierr = SNESDestroy(&snes_power);CHKERRQ(ierr); 673*c4762a1bSJed Brown ierr = SNESDestroy(&snes_water);CHKERRQ(ierr); 674*c4762a1bSJed Brown 675*c4762a1bSJed Brown ierr = DMDestroy(&networkdm);CHKERRQ(ierr); 676*c4762a1bSJed Brown ierr = PetscFinalize(); 677*c4762a1bSJed Brown return ierr; 678*c4762a1bSJed Brown } 679*c4762a1bSJed Brown 680*c4762a1bSJed Brown /*TEST 681*c4762a1bSJed Brown 682*c4762a1bSJed Brown build: 683*c4762a1bSJed Brown requires: !complex double define(PETSC_HAVE_ATTRIBUTEALIGNED) 684*c4762a1bSJed Brown depends: power/PFReadData.c power/pffunctions.c water/waterreaddata.c water/waterfunctions.c 685*c4762a1bSJed Brown 686*c4762a1bSJed Brown test: 687*c4762a1bSJed Brown args: -coupled_snes_converged_reason -options_left no 688*c4762a1bSJed Brown localrunfiles: ex1options power/case9.m water/sample1.inp 689*c4762a1bSJed Brown output_file: output/ex1.out 690*c4762a1bSJed Brown requires: double !complex define(PETSC_HAVE_ATTRIBUTEALIGNED) 691*c4762a1bSJed Brown 692*c4762a1bSJed Brown test: 693*c4762a1bSJed Brown suffix: 2 694*c4762a1bSJed Brown nsize: 3 695*c4762a1bSJed Brown args: -coupled_snes_converged_reason -options_left no 696*c4762a1bSJed Brown localrunfiles: ex1options power/case9.m water/sample1.inp 697*c4762a1bSJed Brown output_file: output/ex1.out 698*c4762a1bSJed Brown requires: double !complex define(PETSC_HAVE_ATTRIBUTEALIGNED) 699*c4762a1bSJed Brown 700*c4762a1bSJed Brown TEST*/ 701