xref: /petsc/src/ts/tutorials/network/pipeInterface.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
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 */
15*9371c9d4SSatish Balay PetscErrorCode PipeCreate(MPI_Comm comm, Pipe *pipe) {
1642dc13f1SHong Zhang   PetscFunctionBegin;
179566063dSJacob Faibussowitsch   PetscCall(PetscNew(pipe));
1842dc13f1SHong Zhang   PetscFunctionReturn(0);
1942dc13f1SHong Zhang }
2042dc13f1SHong Zhang 
2142dc13f1SHong Zhang /*
2242dc13f1SHong Zhang    PipeDestroy - Destroy Pipe object.
2342dc13f1SHong Zhang 
2442dc13f1SHong Zhang    Input Parameters:
2542dc13f1SHong Zhang    pipe - Reference to pipe intended to be destroyed.
2642dc13f1SHong Zhang */
27*9371c9d4SSatish Balay PetscErrorCode PipeDestroy(Pipe *pipe) {
2842dc13f1SHong Zhang   PetscFunctionBegin;
2942dc13f1SHong Zhang   if (!*pipe) PetscFunctionReturn(0);
3042dc13f1SHong Zhang 
319566063dSJacob Faibussowitsch   PetscCall(PipeDestroyJacobian(*pipe));
329566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&(*pipe)->x));
339566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&(*pipe)->da));
3442dc13f1SHong Zhang   PetscFunctionReturn(0);
3542dc13f1SHong Zhang }
3642dc13f1SHong Zhang 
3742dc13f1SHong Zhang /*
3842dc13f1SHong Zhang    PipeSetParameters - Set parameters for Pipe context
3942dc13f1SHong Zhang 
4042dc13f1SHong Zhang    Input Parameter:
4142dc13f1SHong Zhang +  pipe - PIPE object
4242dc13f1SHong Zhang .  length -
4342dc13f1SHong Zhang .  nnodes -
4442dc13f1SHong Zhang .  D -
4542dc13f1SHong Zhang .  a -
4642dc13f1SHong Zhang -  fric -
4742dc13f1SHong Zhang */
48*9371c9d4SSatish Balay PetscErrorCode PipeSetParameters(Pipe pipe, PetscReal length, PetscReal D, PetscReal a, PetscReal fric) {
4942dc13f1SHong Zhang   PetscFunctionBegin;
5042dc13f1SHong Zhang   pipe->length = length;
5142dc13f1SHong Zhang   pipe->D      = D;
5242dc13f1SHong Zhang   pipe->a      = a;
5342dc13f1SHong Zhang   pipe->fric   = fric;
5442dc13f1SHong Zhang   PetscFunctionReturn(0);
5542dc13f1SHong Zhang }
5642dc13f1SHong Zhang 
5742dc13f1SHong Zhang /*
5842dc13f1SHong Zhang     PipeSetUp - Set up pipe based on set parameters.
5942dc13f1SHong Zhang */
60*9371c9d4SSatish Balay PetscErrorCode PipeSetUp(Pipe pipe) {
6142dc13f1SHong Zhang   DMDALocalInfo info;
6242dc13f1SHong Zhang 
6342dc13f1SHong Zhang   PetscFunctionBegin;
649566063dSJacob Faibussowitsch   PetscCall(DMDACreate1d(PETSC_COMM_SELF, DM_BOUNDARY_GHOSTED, pipe->nnodes, 2, 1, NULL, &pipe->da));
659566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(pipe->da));
669566063dSJacob Faibussowitsch   PetscCall(DMSetUp(pipe->da));
679566063dSJacob Faibussowitsch   PetscCall(DMDASetFieldName(pipe->da, 0, "Q"));
689566063dSJacob Faibussowitsch   PetscCall(DMDASetFieldName(pipe->da, 1, "H"));
699566063dSJacob Faibussowitsch   PetscCall(DMDASetUniformCoordinates(pipe->da, 0, pipe->length, 0, 0, 0, 0));
709566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(pipe->da, &(pipe->x)));
7142dc13f1SHong Zhang 
729566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(pipe->da, &info));
7342dc13f1SHong Zhang 
7442dc13f1SHong Zhang   pipe->rad = pipe->D / 2;
7542dc13f1SHong Zhang   pipe->A   = PETSC_PI * pipe->rad * pipe->rad;
7642dc13f1SHong Zhang   pipe->R   = pipe->fric / (2 * pipe->D * pipe->A);
7742dc13f1SHong Zhang   PetscFunctionReturn(0);
7842dc13f1SHong Zhang }
7942dc13f1SHong Zhang 
8042dc13f1SHong Zhang /*
8142dc13f1SHong Zhang     PipeCreateJacobian - Create Jacobian matrix structures for a Pipe.
8242dc13f1SHong Zhang 
8342dc13f1SHong Zhang     Collective on Pipe
8442dc13f1SHong Zhang 
8542dc13f1SHong Zhang     Input Parameter:
8642dc13f1SHong Zhang +   pipe - the Pipe object
8742dc13f1SHong Zhang -   Jin - array of three constructed Jacobian matrices to be reused. Set NULL if it is not available
8842dc13f1SHong Zhang 
8942dc13f1SHong Zhang     Output Parameter:
9042dc13f1SHong Zhang .   J  - array of three empty Jacobian matrices
9142dc13f1SHong Zhang 
9242dc13f1SHong Zhang     Level: beginner
9342dc13f1SHong Zhang */
94*9371c9d4SSatish Balay PetscErrorCode PipeCreateJacobian(Pipe pipe, Mat *Jin, Mat *J[]) {
9542dc13f1SHong Zhang   Mat         *Jpipe;
9642dc13f1SHong Zhang   PetscInt     M, rows[2], cols[2], *nz;
9742dc13f1SHong Zhang   PetscScalar *aa;
9842dc13f1SHong Zhang 
9942dc13f1SHong Zhang   PetscFunctionBegin;
10042dc13f1SHong Zhang   if (Jin) {
10142dc13f1SHong Zhang     *J             = Jin;
10242dc13f1SHong Zhang     pipe->jacobian = Jin;
1039566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)(Jin[0])));
10442dc13f1SHong Zhang     PetscFunctionReturn(0);
10542dc13f1SHong Zhang   }
1069566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(3, &Jpipe));
10742dc13f1SHong Zhang 
10842dc13f1SHong Zhang   /* Jacobian for this pipe */
1099566063dSJacob Faibussowitsch   PetscCall(DMSetMatrixStructureOnly(pipe->da, PETSC_TRUE));
1109566063dSJacob Faibussowitsch   PetscCall(DMCreateMatrix(pipe->da, &Jpipe[0]));
1119566063dSJacob Faibussowitsch   PetscCall(DMSetMatrixStructureOnly(pipe->da, PETSC_FALSE));
11242dc13f1SHong Zhang 
11342dc13f1SHong Zhang   /* Jacobian for upstream vertex */
1149566063dSJacob Faibussowitsch   PetscCall(MatGetSize(Jpipe[0], &M, NULL));
1159566063dSJacob Faibussowitsch   PetscCall(PetscCalloc2(M, &nz, 4, &aa));
11642dc13f1SHong Zhang 
1179566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_SELF, &Jpipe[1]));
1189566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(Jpipe[1], PETSC_DECIDE, PETSC_DECIDE, M, 2));
1199566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(Jpipe[1]));
1209566063dSJacob Faibussowitsch   PetscCall(MatSetOption(Jpipe[1], MAT_STRUCTURE_ONLY, PETSC_TRUE));
121*9371c9d4SSatish Balay   nz[0]   = 2;
122*9371c9d4SSatish Balay   nz[1]   = 2;
123*9371c9d4SSatish Balay   rows[0] = 0;
124*9371c9d4SSatish Balay   rows[1] = 1;
125*9371c9d4SSatish Balay   cols[0] = 0;
126*9371c9d4SSatish Balay   cols[1] = 1;
1279566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(Jpipe[1], 0, nz));
1289566063dSJacob Faibussowitsch   PetscCall(MatSetValues(Jpipe[1], 2, rows, 2, cols, aa, INSERT_VALUES));
1299566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(Jpipe[1], MAT_FINAL_ASSEMBLY));
1309566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(Jpipe[1], MAT_FINAL_ASSEMBLY));
13142dc13f1SHong Zhang 
13242dc13f1SHong Zhang   /* Jacobian for downstream vertex */
1339566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_SELF, &Jpipe[2]));
1349566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(Jpipe[2], PETSC_DECIDE, PETSC_DECIDE, M, 2));
1359566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(Jpipe[2]));
1369566063dSJacob Faibussowitsch   PetscCall(MatSetOption(Jpipe[2], MAT_STRUCTURE_ONLY, PETSC_TRUE));
137*9371c9d4SSatish Balay   nz[0]     = 0;
138*9371c9d4SSatish Balay   nz[1]     = 0;
139*9371c9d4SSatish Balay   nz[M - 2] = 2;
140*9371c9d4SSatish Balay   nz[M - 1] = 2;
141*9371c9d4SSatish Balay   rows[0]   = M - 2;
142*9371c9d4SSatish Balay   rows[1]   = M - 1;
1439566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(Jpipe[2], 0, nz));
1449566063dSJacob Faibussowitsch   PetscCall(MatSetValues(Jpipe[2], 2, rows, 2, cols, aa, INSERT_VALUES));
1459566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(Jpipe[2], MAT_FINAL_ASSEMBLY));
1469566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(Jpipe[2], MAT_FINAL_ASSEMBLY));
14742dc13f1SHong Zhang 
1489566063dSJacob Faibussowitsch   PetscCall(PetscFree2(nz, aa));
14942dc13f1SHong Zhang 
15042dc13f1SHong Zhang   *J             = Jpipe;
15142dc13f1SHong Zhang   pipe->jacobian = Jpipe;
15242dc13f1SHong Zhang   PetscFunctionReturn(0);
15342dc13f1SHong Zhang }
15442dc13f1SHong Zhang 
155*9371c9d4SSatish Balay PetscErrorCode PipeDestroyJacobian(Pipe pipe) {
15642dc13f1SHong Zhang   Mat     *Jpipe = pipe->jacobian;
15742dc13f1SHong Zhang   PetscInt i;
15842dc13f1SHong Zhang 
15942dc13f1SHong Zhang   PetscFunctionBegin;
16042dc13f1SHong Zhang   if (Jpipe) {
161*9371c9d4SSatish Balay     for (i = 0; i < 3; i++) { PetscCall(MatDestroy(&Jpipe[i])); }
16242dc13f1SHong Zhang   }
1639566063dSJacob Faibussowitsch   PetscCall(PetscFree(Jpipe));
16442dc13f1SHong Zhang   PetscFunctionReturn(0);
16542dc13f1SHong Zhang }
16642dc13f1SHong Zhang 
16742dc13f1SHong Zhang /*
16842dc13f1SHong Zhang     JunctionCreateJacobian - Create Jacobian matrices for a vertex.
16942dc13f1SHong Zhang 
17042dc13f1SHong Zhang     Collective on Pipe
17142dc13f1SHong Zhang 
17242dc13f1SHong Zhang     Input Parameter:
17342dc13f1SHong Zhang +   dm - the DMNetwork object
17442dc13f1SHong Zhang .   v - vertex point
17542dc13f1SHong Zhang -   Jin - Jacobian patterns created by JunctionCreateJacobianSample() for reuse
17642dc13f1SHong Zhang 
17742dc13f1SHong Zhang     Output Parameter:
17842dc13f1SHong Zhang .   J  - array of Jacobian matrices (see dmnetworkimpl.h)
17942dc13f1SHong Zhang 
18042dc13f1SHong Zhang     Level: beginner
18142dc13f1SHong Zhang */
182*9371c9d4SSatish Balay PetscErrorCode JunctionCreateJacobian(DM dm, PetscInt v, Mat *Jin, Mat *J[]) {
18342dc13f1SHong Zhang   Mat            *Jv;
18442dc13f1SHong Zhang   PetscInt        nedges, e, i, M, N, *rows, *cols;
18542dc13f1SHong Zhang   PetscBool       isSelf;
18642dc13f1SHong Zhang   const PetscInt *edges, *cone;
18742dc13f1SHong Zhang   PetscScalar    *zeros;
18842dc13f1SHong Zhang 
18942dc13f1SHong Zhang   PetscFunctionBegin;
190a5b23f4aSJose E. Roman   /* Get array size of Jv */
1919566063dSJacob Faibussowitsch   PetscCall(DMNetworkGetSupportingEdges(dm, v, &nedges, &edges));
19263a3b9bcSJacob Faibussowitsch   PetscCheck(nedges > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "%" PetscInt_FMT " vertex, nedges %" PetscInt_FMT, v, nedges);
19342dc13f1SHong Zhang 
19442dc13f1SHong Zhang   /* two Jacobians for each connected edge: J(v,e) and J(v,vc); adding J(v,v), total 2*nedges+1 Jacobians */
1959566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(2 * nedges + 1, &Jv));
19642dc13f1SHong Zhang 
19742dc13f1SHong Zhang   /* Create dense zero block for this vertex: J[0] = Jacobian(v,v) */
1989566063dSJacob Faibussowitsch   PetscCall(DMNetworkGetComponent(dm, v, -1, NULL, NULL, &M));
19963a3b9bcSJacob Faibussowitsch   PetscCheck(M == 2, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "%" PetscInt_FMT " != 2", M);
2009566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(M, &rows, M, &cols, M * M, &zeros));
2019566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(zeros, M * M));
20242dc13f1SHong Zhang   for (i = 0; i < M; i++) rows[i] = i;
20342dc13f1SHong Zhang 
20442dc13f1SHong Zhang   for (e = 0; e < nedges; e++) {
20542dc13f1SHong Zhang     /* create Jv[2*e+1] = Jacobian(v,e), e: supporting edge */
2069566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetConnectedVertices(dm, edges[e], &cone));
20742dc13f1SHong Zhang     isSelf = (v == cone[0]) ? PETSC_TRUE : PETSC_FALSE;
20842dc13f1SHong Zhang 
20942dc13f1SHong Zhang     if (Jin) {
21042dc13f1SHong Zhang       if (isSelf) {
21142dc13f1SHong Zhang         Jv[2 * e + 1] = Jin[0];
21242dc13f1SHong Zhang       } else {
21342dc13f1SHong Zhang         Jv[2 * e + 1] = Jin[1];
21442dc13f1SHong Zhang       }
21542dc13f1SHong Zhang       Jv[2 * e + 2] = Jin[2];
2169566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)(Jv[2 * e + 1])));
2179566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)(Jv[2 * e + 2])));
21842dc13f1SHong Zhang     } else {
21942dc13f1SHong Zhang       /* create J(v,e) */
2209566063dSJacob Faibussowitsch       PetscCall(MatCreate(PETSC_COMM_SELF, &Jv[2 * e + 1]));
2219566063dSJacob Faibussowitsch       PetscCall(DMNetworkGetComponent(dm, edges[e], -1, NULL, NULL, &N));
2229566063dSJacob Faibussowitsch       PetscCall(MatSetSizes(Jv[2 * e + 1], PETSC_DECIDE, PETSC_DECIDE, M, N));
2239566063dSJacob Faibussowitsch       PetscCall(MatSetFromOptions(Jv[2 * e + 1]));
2249566063dSJacob Faibussowitsch       PetscCall(MatSetOption(Jv[2 * e + 1], MAT_STRUCTURE_ONLY, PETSC_TRUE));
2259566063dSJacob Faibussowitsch       PetscCall(MatSeqAIJSetPreallocation(Jv[2 * e + 1], 2, NULL));
22642dc13f1SHong Zhang       if (N) {
22742dc13f1SHong Zhang         if (isSelf) { /* coupling at upstream */
22842dc13f1SHong Zhang           for (i = 0; i < 2; i++) cols[i] = i;
22942dc13f1SHong Zhang         } else { /* coupling at downstream */
230*9371c9d4SSatish Balay           cols[0] = N - 2;
231*9371c9d4SSatish Balay           cols[1] = N - 1;
23242dc13f1SHong Zhang         }
2339566063dSJacob Faibussowitsch         PetscCall(MatSetValues(Jv[2 * e + 1], 2, rows, 2, cols, zeros, INSERT_VALUES));
23442dc13f1SHong Zhang       }
2359566063dSJacob Faibussowitsch       PetscCall(MatAssemblyBegin(Jv[2 * e + 1], MAT_FINAL_ASSEMBLY));
2369566063dSJacob Faibussowitsch       PetscCall(MatAssemblyEnd(Jv[2 * e + 1], MAT_FINAL_ASSEMBLY));
23742dc13f1SHong Zhang 
23842dc13f1SHong Zhang       /* create Jv[2*e+2] = Jacobian(v,vc), vc: connected vertex.
23942dc13f1SHong Zhang        In WashNetwork, v and vc are not connected, thus Jacobian(v,vc) is empty */
2409566063dSJacob Faibussowitsch       PetscCall(MatCreate(PETSC_COMM_SELF, &Jv[2 * e + 2]));
2419566063dSJacob Faibussowitsch       PetscCall(MatSetSizes(Jv[2 * e + 2], PETSC_DECIDE, PETSC_DECIDE, M, M)); /* empty matrix, sizes can be arbitrary */
2429566063dSJacob Faibussowitsch       PetscCall(MatSetFromOptions(Jv[2 * e + 2]));
2439566063dSJacob Faibussowitsch       PetscCall(MatSetOption(Jv[2 * e + 2], MAT_STRUCTURE_ONLY, PETSC_TRUE));
2449566063dSJacob Faibussowitsch       PetscCall(MatSeqAIJSetPreallocation(Jv[2 * e + 2], 1, NULL));
2459566063dSJacob Faibussowitsch       PetscCall(MatAssemblyBegin(Jv[2 * e + 2], MAT_FINAL_ASSEMBLY));
2469566063dSJacob Faibussowitsch       PetscCall(MatAssemblyEnd(Jv[2 * e + 2], MAT_FINAL_ASSEMBLY));
24742dc13f1SHong Zhang     }
24842dc13f1SHong Zhang   }
2499566063dSJacob Faibussowitsch   PetscCall(PetscFree3(rows, cols, zeros));
25042dc13f1SHong Zhang 
25142dc13f1SHong Zhang   *J = Jv;
25242dc13f1SHong Zhang   PetscFunctionReturn(0);
25342dc13f1SHong Zhang }
25442dc13f1SHong Zhang 
255*9371c9d4SSatish Balay PetscErrorCode JunctionDestroyJacobian(DM dm, PetscInt v, Junction junc) {
25642dc13f1SHong Zhang   Mat            *Jv = junc->jacobian;
25742dc13f1SHong Zhang   const PetscInt *edges;
25842dc13f1SHong Zhang   PetscInt        nedges, e;
25942dc13f1SHong Zhang 
26042dc13f1SHong Zhang   PetscFunctionBegin;
26142dc13f1SHong Zhang   if (!Jv) PetscFunctionReturn(0);
26242dc13f1SHong Zhang 
2639566063dSJacob Faibussowitsch   PetscCall(DMNetworkGetSupportingEdges(dm, v, &nedges, &edges));
26442dc13f1SHong Zhang   for (e = 0; e < nedges; e++) {
2659566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&Jv[2 * e + 1]));
2669566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&Jv[2 * e + 2]));
26742dc13f1SHong Zhang   }
2689566063dSJacob Faibussowitsch   PetscCall(PetscFree(Jv));
26942dc13f1SHong Zhang   PetscFunctionReturn(0);
27042dc13f1SHong Zhang }
271