xref: /petsc/src/ts/tutorials/network/pipeInterface.c (revision 9566063d113dddea24716c546802770db7481bc0)
142dc13f1SHong Zhang #include "wash.h"
242dc13f1SHong Zhang 
342dc13f1SHong Zhang /* Subroutines for Pipe                                  */
442dc13f1SHong Zhang /* -------------------------------------------------------*/
542dc13f1SHong Zhang 
642dc13f1SHong Zhang /*
742dc13f1SHong Zhang    PipeCreate - Create Pipe object.
842dc13f1SHong Zhang 
942dc13f1SHong Zhang    Input Parameters:
1042dc13f1SHong Zhang    comm - MPI communicator
1142dc13f1SHong Zhang 
1242dc13f1SHong Zhang    Output Parameter:
1342dc13f1SHong Zhang .  pipe - location to put the PIPE context
1442dc13f1SHong Zhang */
1542dc13f1SHong Zhang PetscErrorCode PipeCreate(MPI_Comm comm,Pipe *pipe)
1642dc13f1SHong Zhang {
1742dc13f1SHong Zhang   PetscFunctionBegin;
18*9566063dSJacob Faibussowitsch   PetscCall(PetscNew(pipe));
1942dc13f1SHong Zhang   PetscFunctionReturn(0);
2042dc13f1SHong Zhang }
2142dc13f1SHong Zhang 
2242dc13f1SHong Zhang /*
2342dc13f1SHong Zhang    PipeDestroy - Destroy Pipe object.
2442dc13f1SHong Zhang 
2542dc13f1SHong Zhang    Input Parameters:
2642dc13f1SHong Zhang    pipe - Reference to pipe intended to be destroyed.
2742dc13f1SHong Zhang */
2842dc13f1SHong Zhang PetscErrorCode PipeDestroy(Pipe *pipe)
2942dc13f1SHong Zhang {
3042dc13f1SHong Zhang   PetscFunctionBegin;
3142dc13f1SHong Zhang   if (!*pipe) PetscFunctionReturn(0);
3242dc13f1SHong Zhang 
33*9566063dSJacob Faibussowitsch   PetscCall(PipeDestroyJacobian(*pipe));
34*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&(*pipe)->x));
35*9566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&(*pipe)->da));
3642dc13f1SHong Zhang   PetscFunctionReturn(0);
3742dc13f1SHong Zhang }
3842dc13f1SHong Zhang 
3942dc13f1SHong Zhang /*
4042dc13f1SHong Zhang    PipeSetParameters - Set parameters for Pipe context
4142dc13f1SHong Zhang 
4242dc13f1SHong Zhang    Input Parameter:
4342dc13f1SHong Zhang +  pipe - PIPE object
4442dc13f1SHong Zhang .  length -
4542dc13f1SHong Zhang .  nnodes -
4642dc13f1SHong Zhang .  D -
4742dc13f1SHong Zhang .  a -
4842dc13f1SHong Zhang -  fric -
4942dc13f1SHong Zhang */
5042dc13f1SHong Zhang PetscErrorCode PipeSetParameters(Pipe pipe,PetscReal length,PetscReal D,PetscReal a,PetscReal fric)
5142dc13f1SHong Zhang {
5242dc13f1SHong Zhang   PetscFunctionBegin;
5342dc13f1SHong Zhang   pipe->length = length;
5442dc13f1SHong Zhang   pipe->D      = D;
5542dc13f1SHong Zhang   pipe->a      = a;
5642dc13f1SHong Zhang   pipe->fric   = fric;
5742dc13f1SHong Zhang   PetscFunctionReturn(0);
5842dc13f1SHong Zhang }
5942dc13f1SHong Zhang 
6042dc13f1SHong Zhang /*
6142dc13f1SHong Zhang     PipeSetUp - Set up pipe based on set parameters.
6242dc13f1SHong Zhang */
6342dc13f1SHong Zhang PetscErrorCode PipeSetUp(Pipe pipe)
6442dc13f1SHong Zhang {
6542dc13f1SHong Zhang   DMDALocalInfo  info;
6642dc13f1SHong Zhang 
6742dc13f1SHong Zhang   PetscFunctionBegin;
68*9566063dSJacob Faibussowitsch   PetscCall(DMDACreate1d(PETSC_COMM_SELF, DM_BOUNDARY_GHOSTED, pipe->nnodes, 2, 1, NULL, &pipe->da));
69*9566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(pipe->da));
70*9566063dSJacob Faibussowitsch   PetscCall(DMSetUp(pipe->da));
71*9566063dSJacob Faibussowitsch   PetscCall(DMDASetFieldName(pipe->da, 0, "Q"));
72*9566063dSJacob Faibussowitsch   PetscCall(DMDASetFieldName(pipe->da, 1, "H"));
73*9566063dSJacob Faibussowitsch   PetscCall(DMDASetUniformCoordinates(pipe->da, 0, pipe->length, 0, 0, 0, 0));
74*9566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(pipe->da, &(pipe->x)));
7542dc13f1SHong Zhang 
76*9566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(pipe->da, &info));
7742dc13f1SHong Zhang 
7842dc13f1SHong Zhang   pipe->rad = pipe->D / 2;
7942dc13f1SHong Zhang   pipe->A   = PETSC_PI*pipe->rad*pipe->rad;
8042dc13f1SHong Zhang   pipe->R   = pipe->fric / (2*pipe->D*pipe->A);
8142dc13f1SHong Zhang   PetscFunctionReturn(0);
8242dc13f1SHong Zhang }
8342dc13f1SHong Zhang 
8442dc13f1SHong Zhang /*
8542dc13f1SHong Zhang     PipeCreateJacobian - Create Jacobian matrix structures for a Pipe.
8642dc13f1SHong Zhang 
8742dc13f1SHong Zhang     Collective on Pipe
8842dc13f1SHong Zhang 
8942dc13f1SHong Zhang     Input Parameter:
9042dc13f1SHong Zhang +   pipe - the Pipe object
9142dc13f1SHong Zhang -   Jin - array of three constructed Jacobian matrices to be reused. Set NULL if it is not available
9242dc13f1SHong Zhang 
9342dc13f1SHong Zhang     Output Parameter:
9442dc13f1SHong Zhang .   J  - array of three empty Jacobian matrices
9542dc13f1SHong Zhang 
9642dc13f1SHong Zhang     Level: beginner
9742dc13f1SHong Zhang */
9842dc13f1SHong Zhang PetscErrorCode PipeCreateJacobian(Pipe pipe,Mat *Jin,Mat *J[])
9942dc13f1SHong Zhang {
10042dc13f1SHong Zhang   Mat            *Jpipe;
10142dc13f1SHong Zhang   PetscInt       M,rows[2],cols[2],*nz;
10242dc13f1SHong Zhang   PetscScalar    *aa;
10342dc13f1SHong Zhang 
10442dc13f1SHong Zhang   PetscFunctionBegin;
10542dc13f1SHong Zhang   if (Jin) {
10642dc13f1SHong Zhang     *J = Jin;
10742dc13f1SHong Zhang     pipe->jacobian = Jin;
108*9566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)(Jin[0])));
10942dc13f1SHong Zhang     PetscFunctionReturn(0);
11042dc13f1SHong Zhang   }
111*9566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(3,&Jpipe));
11242dc13f1SHong Zhang 
11342dc13f1SHong Zhang   /* Jacobian for this pipe */
114*9566063dSJacob Faibussowitsch   PetscCall(DMSetMatrixStructureOnly(pipe->da,PETSC_TRUE));
115*9566063dSJacob Faibussowitsch   PetscCall(DMCreateMatrix(pipe->da,&Jpipe[0]));
116*9566063dSJacob Faibussowitsch   PetscCall(DMSetMatrixStructureOnly(pipe->da,PETSC_FALSE));
11742dc13f1SHong Zhang 
11842dc13f1SHong Zhang   /* Jacobian for upstream vertex */
119*9566063dSJacob Faibussowitsch   PetscCall(MatGetSize(Jpipe[0],&M,NULL));
120*9566063dSJacob Faibussowitsch   PetscCall(PetscCalloc2(M,&nz,4,&aa));
12142dc13f1SHong Zhang 
122*9566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_SELF,&Jpipe[1]));
123*9566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(Jpipe[1],PETSC_DECIDE,PETSC_DECIDE,M,2));
124*9566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(Jpipe[1]));
125*9566063dSJacob Faibussowitsch   PetscCall(MatSetOption(Jpipe[1],MAT_STRUCTURE_ONLY,PETSC_TRUE));
12642dc13f1SHong Zhang   nz[0] = 2; nz[1] = 2;
12742dc13f1SHong Zhang   rows[0] = 0; rows[1] = 1;
12842dc13f1SHong Zhang   cols[0] = 0; cols[1] = 1;
129*9566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(Jpipe[1],0,nz));
130*9566063dSJacob Faibussowitsch   PetscCall(MatSetValues(Jpipe[1],2,rows,2,cols,aa,INSERT_VALUES));
131*9566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(Jpipe[1],MAT_FINAL_ASSEMBLY));
132*9566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(Jpipe[1],MAT_FINAL_ASSEMBLY));
13342dc13f1SHong Zhang 
13442dc13f1SHong Zhang   /* Jacobian for downstream vertex */
135*9566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_SELF,&Jpipe[2]));
136*9566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(Jpipe[2],PETSC_DECIDE,PETSC_DECIDE,M,2));
137*9566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(Jpipe[2]));
138*9566063dSJacob Faibussowitsch   PetscCall(MatSetOption(Jpipe[2],MAT_STRUCTURE_ONLY,PETSC_TRUE));
13942dc13f1SHong Zhang   nz[0] = 0; nz[1] = 0; nz[M-2] = 2; nz[M-1] = 2;
14042dc13f1SHong Zhang   rows[0] = M - 2; rows[1] = M - 1;
141*9566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(Jpipe[2],0,nz));
142*9566063dSJacob Faibussowitsch   PetscCall(MatSetValues(Jpipe[2],2,rows,2,cols,aa,INSERT_VALUES));
143*9566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(Jpipe[2],MAT_FINAL_ASSEMBLY));
144*9566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(Jpipe[2],MAT_FINAL_ASSEMBLY));
14542dc13f1SHong Zhang 
146*9566063dSJacob Faibussowitsch   PetscCall(PetscFree2(nz,aa));
14742dc13f1SHong Zhang 
14842dc13f1SHong Zhang   *J = Jpipe;
14942dc13f1SHong Zhang   pipe->jacobian = Jpipe;
15042dc13f1SHong Zhang   PetscFunctionReturn(0);
15142dc13f1SHong Zhang }
15242dc13f1SHong Zhang 
15342dc13f1SHong Zhang PetscErrorCode PipeDestroyJacobian(Pipe pipe)
15442dc13f1SHong Zhang {
15542dc13f1SHong Zhang   Mat            *Jpipe = pipe->jacobian;
15642dc13f1SHong Zhang   PetscInt       i;
15742dc13f1SHong Zhang 
15842dc13f1SHong Zhang   PetscFunctionBegin;
15942dc13f1SHong Zhang   if (Jpipe) {
16042dc13f1SHong Zhang     for (i=0; i<3; i++) {
161*9566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&Jpipe[i]));
16242dc13f1SHong Zhang     }
16342dc13f1SHong Zhang   }
164*9566063dSJacob Faibussowitsch   PetscCall(PetscFree(Jpipe));
16542dc13f1SHong Zhang   PetscFunctionReturn(0);
16642dc13f1SHong Zhang }
16742dc13f1SHong Zhang 
16842dc13f1SHong Zhang /*
16942dc13f1SHong Zhang     JunctionCreateJacobian - Create Jacobian matrices for a vertex.
17042dc13f1SHong Zhang 
17142dc13f1SHong Zhang     Collective on Pipe
17242dc13f1SHong Zhang 
17342dc13f1SHong Zhang     Input Parameter:
17442dc13f1SHong Zhang +   dm - the DMNetwork object
17542dc13f1SHong Zhang .   v - vertex point
17642dc13f1SHong Zhang -   Jin - Jacobian patterns created by JunctionCreateJacobianSample() for reuse
17742dc13f1SHong Zhang 
17842dc13f1SHong Zhang     Output Parameter:
17942dc13f1SHong Zhang .   J  - array of Jacobian matrices (see dmnetworkimpl.h)
18042dc13f1SHong Zhang 
18142dc13f1SHong Zhang     Level: beginner
18242dc13f1SHong Zhang */
18342dc13f1SHong Zhang PetscErrorCode JunctionCreateJacobian(DM dm,PetscInt v,Mat *Jin,Mat *J[])
18442dc13f1SHong Zhang {
18542dc13f1SHong Zhang   Mat            *Jv;
18642dc13f1SHong Zhang   PetscInt       nedges,e,i,M,N,*rows,*cols;
18742dc13f1SHong Zhang   PetscBool      isSelf;
18842dc13f1SHong Zhang   const PetscInt *edges,*cone;
18942dc13f1SHong Zhang   PetscScalar    *zeros;
19042dc13f1SHong Zhang 
19142dc13f1SHong Zhang   PetscFunctionBegin;
192a5b23f4aSJose E. Roman   /* Get array size of Jv */
193*9566063dSJacob Faibussowitsch   PetscCall(DMNetworkGetSupportingEdges(dm,v,&nedges,&edges));
1942c71b3e2SJacob Faibussowitsch   PetscCheckFalse(nedges <= 0,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"%d vertex, nedges %d",v,nedges);
19542dc13f1SHong Zhang 
19642dc13f1SHong Zhang   /* two Jacobians for each connected edge: J(v,e) and J(v,vc); adding J(v,v), total 2*nedges+1 Jacobians */
197*9566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(2*nedges+1,&Jv));
19842dc13f1SHong Zhang 
19942dc13f1SHong Zhang   /* Create dense zero block for this vertex: J[0] = Jacobian(v,v) */
200*9566063dSJacob Faibussowitsch   PetscCall(DMNetworkGetComponent(dm,v,-1,NULL,NULL,&M));
2012c71b3e2SJacob Faibussowitsch   PetscCheckFalse(M !=2,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"M != 2",M);
202*9566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(M,&rows,M,&cols,M*M,&zeros));
203*9566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(zeros,M*M));
20442dc13f1SHong Zhang   for (i=0; i<M; i++) rows[i] = i;
20542dc13f1SHong Zhang 
20642dc13f1SHong Zhang   for (e=0; e<nedges; e++) {
20742dc13f1SHong Zhang     /* create Jv[2*e+1] = Jacobian(v,e), e: supporting edge */
208*9566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetConnectedVertices(dm,edges[e],&cone));
20942dc13f1SHong Zhang     isSelf = (v == cone[0]) ? PETSC_TRUE:PETSC_FALSE;
21042dc13f1SHong Zhang 
21142dc13f1SHong Zhang     if (Jin) {
21242dc13f1SHong Zhang       if (isSelf) {
21342dc13f1SHong Zhang         Jv[2*e+1] = Jin[0];
21442dc13f1SHong Zhang       } else {
21542dc13f1SHong Zhang         Jv[2*e+1] = Jin[1];
21642dc13f1SHong Zhang       }
21742dc13f1SHong Zhang       Jv[2*e+2] = Jin[2];
218*9566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)(Jv[2*e+1])));
219*9566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)(Jv[2*e+2])));
22042dc13f1SHong Zhang     } else {
22142dc13f1SHong Zhang       /* create J(v,e) */
222*9566063dSJacob Faibussowitsch       PetscCall(MatCreate(PETSC_COMM_SELF,&Jv[2*e+1]));
223*9566063dSJacob Faibussowitsch       PetscCall(DMNetworkGetComponent(dm,edges[e],-1,NULL,NULL,&N));
224*9566063dSJacob Faibussowitsch       PetscCall(MatSetSizes(Jv[2*e+1],PETSC_DECIDE,PETSC_DECIDE,M,N));
225*9566063dSJacob Faibussowitsch       PetscCall(MatSetFromOptions(Jv[2*e+1]));
226*9566063dSJacob Faibussowitsch       PetscCall(MatSetOption(Jv[2*e+1],MAT_STRUCTURE_ONLY,PETSC_TRUE));
227*9566063dSJacob Faibussowitsch       PetscCall(MatSeqAIJSetPreallocation(Jv[2*e+1],2,NULL));
22842dc13f1SHong Zhang       if (N) {
22942dc13f1SHong Zhang         if (isSelf) { /* coupling at upstream */
23042dc13f1SHong Zhang           for (i=0; i<2; i++) cols[i] = i;
23142dc13f1SHong Zhang         } else { /* coupling at downstream */
23242dc13f1SHong Zhang           cols[0] = N-2; cols[1] = N-1;
23342dc13f1SHong Zhang         }
234*9566063dSJacob Faibussowitsch         PetscCall(MatSetValues(Jv[2*e+1],2,rows,2,cols,zeros,INSERT_VALUES));
23542dc13f1SHong Zhang       }
236*9566063dSJacob Faibussowitsch       PetscCall(MatAssemblyBegin(Jv[2*e+1],MAT_FINAL_ASSEMBLY));
237*9566063dSJacob Faibussowitsch       PetscCall(MatAssemblyEnd(Jv[2*e+1],MAT_FINAL_ASSEMBLY));
23842dc13f1SHong Zhang 
23942dc13f1SHong Zhang       /* create Jv[2*e+2] = Jacobian(v,vc), vc: connected vertex.
24042dc13f1SHong Zhang        In WashNetwork, v and vc are not connected, thus Jacobian(v,vc) is empty */
241*9566063dSJacob Faibussowitsch       PetscCall(MatCreate(PETSC_COMM_SELF,&Jv[2*e+2]));
242*9566063dSJacob Faibussowitsch       PetscCall(MatSetSizes(Jv[2*e+2],PETSC_DECIDE,PETSC_DECIDE,M,M)); /* empty matrix, sizes can be arbitrary */
243*9566063dSJacob Faibussowitsch       PetscCall(MatSetFromOptions(Jv[2*e+2]));
244*9566063dSJacob Faibussowitsch       PetscCall(MatSetOption(Jv[2*e+2],MAT_STRUCTURE_ONLY,PETSC_TRUE));
245*9566063dSJacob Faibussowitsch       PetscCall(MatSeqAIJSetPreallocation(Jv[2*e+2],1,NULL));
246*9566063dSJacob Faibussowitsch       PetscCall(MatAssemblyBegin(Jv[2*e+2],MAT_FINAL_ASSEMBLY));
247*9566063dSJacob Faibussowitsch       PetscCall(MatAssemblyEnd(Jv[2*e+2],MAT_FINAL_ASSEMBLY));
24842dc13f1SHong Zhang     }
24942dc13f1SHong Zhang   }
250*9566063dSJacob Faibussowitsch   PetscCall(PetscFree3(rows,cols,zeros));
25142dc13f1SHong Zhang 
25242dc13f1SHong Zhang   *J = Jv;
25342dc13f1SHong Zhang   PetscFunctionReturn(0);
25442dc13f1SHong Zhang }
25542dc13f1SHong Zhang 
25642dc13f1SHong Zhang PetscErrorCode JunctionDestroyJacobian(DM dm,PetscInt v,Junction junc)
25742dc13f1SHong Zhang {
25842dc13f1SHong Zhang   Mat            *Jv=junc->jacobian;
25942dc13f1SHong Zhang   const PetscInt *edges;
26042dc13f1SHong Zhang   PetscInt       nedges,e;
26142dc13f1SHong Zhang 
26242dc13f1SHong Zhang   PetscFunctionBegin;
26342dc13f1SHong Zhang   if (!Jv) PetscFunctionReturn(0);
26442dc13f1SHong Zhang 
265*9566063dSJacob Faibussowitsch   PetscCall(DMNetworkGetSupportingEdges(dm,v,&nedges,&edges));
26642dc13f1SHong Zhang   for (e=0; e<nedges; e++) {
267*9566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&Jv[2*e+1]));
268*9566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&Jv[2*e+2]));
26942dc13f1SHong Zhang   }
270*9566063dSJacob Faibussowitsch   PetscCall(PetscFree(Jv));
27142dc13f1SHong Zhang   PetscFunctionReturn(0);
27242dc13f1SHong Zhang }
273