1c4762a1bSJed Brown static char help[] = "Multiphase flow in a porous medium in 1d.\n\n";
2c4762a1bSJed Brown #include <petscdm.h>
3c4762a1bSJed Brown #include <petscdmda.h>
4c4762a1bSJed Brown #include <petscsnes.h>
5c4762a1bSJed Brown
6c4762a1bSJed Brown typedef struct {
7c4762a1bSJed Brown DM cda;
8c4762a1bSJed Brown Vec uold;
9c4762a1bSJed Brown Vec Kappa;
10c4762a1bSJed Brown PetscReal phi;
11c4762a1bSJed Brown PetscReal kappaWet;
12c4762a1bSJed Brown PetscReal kappaNoWet;
13c4762a1bSJed Brown PetscReal dt;
14c4762a1bSJed Brown /* Boundary conditions */
15c4762a1bSJed Brown PetscReal sl, vl, pl;
16c4762a1bSJed Brown } AppCtx;
17c4762a1bSJed Brown
18c4762a1bSJed Brown typedef struct {
19c4762a1bSJed Brown PetscScalar s; /* The saturation on each cell */
20c4762a1bSJed Brown PetscScalar v; /* The velocity on each face */
21c4762a1bSJed Brown PetscScalar p; /* The pressure on each cell */
22c4762a1bSJed Brown } Field;
23c4762a1bSJed Brown
24c4762a1bSJed Brown /*
25c4762a1bSJed Brown FormPermeability - Forms permeability field.
26c4762a1bSJed Brown
27c4762a1bSJed Brown Input Parameters:
28c4762a1bSJed Brown user - user-defined application context
29c4762a1bSJed Brown
30c4762a1bSJed Brown Output Parameter:
31c4762a1bSJed Brown Kappa - vector
32c4762a1bSJed Brown */
FormPermeability(DM da,Vec Kappa,AppCtx * user)33d71ae5a4SJacob Faibussowitsch PetscErrorCode FormPermeability(DM da, Vec Kappa, AppCtx *user)
34d71ae5a4SJacob Faibussowitsch {
35c4762a1bSJed Brown DM cda;
36c4762a1bSJed Brown Vec c;
37c4762a1bSJed Brown PetscScalar *K;
38c4762a1bSJed Brown PetscScalar *coords;
39c4762a1bSJed Brown PetscInt xs, xm, i;
40c4762a1bSJed Brown
41c4762a1bSJed Brown PetscFunctionBeginUser;
429566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(da, &cda));
439566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(da, &c));
449566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, NULL, NULL, &xm, NULL, NULL));
459566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, Kappa, &K));
469566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(cda, c, &coords));
47c4762a1bSJed Brown for (i = xs; i < xs + xm; ++i) {
48c4762a1bSJed Brown #if 1
49c4762a1bSJed Brown K[i] = 1.0;
50c4762a1bSJed Brown #else
51c4762a1bSJed Brown /* Notch */
52c4762a1bSJed Brown if (i == (xs + xm) / 2) K[i] = 0.00000001;
53c4762a1bSJed Brown else K[i] = 1.0;
54c4762a1bSJed Brown #endif
55c4762a1bSJed Brown }
569566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, Kappa, &K));
579566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(cda, c, &coords));
58*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
59c4762a1bSJed Brown }
60c4762a1bSJed Brown
61c4762a1bSJed Brown /*
62c4762a1bSJed Brown FormFunctionLocal - Evaluates nonlinear residual, F(x) on local process patch
63c4762a1bSJed Brown */
FormFunctionLocal(DMDALocalInfo * info,Field * u,Field * f,AppCtx * user)64d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunctionLocal(DMDALocalInfo *info, Field *u, Field *f, AppCtx *user)
65d71ae5a4SJacob Faibussowitsch {
66c4762a1bSJed Brown Vec L;
67c4762a1bSJed Brown PetscReal phi = user->phi;
68c4762a1bSJed Brown PetscReal dt = user->dt;
69c4762a1bSJed Brown PetscReal dx = 1.0 / (PetscReal)(info->mx - 1);
70c4762a1bSJed Brown PetscReal alpha = 2.0;
71c4762a1bSJed Brown PetscReal beta = 2.0;
72c4762a1bSJed Brown PetscReal kappaWet = user->kappaWet;
73c4762a1bSJed Brown PetscReal kappaNoWet = user->kappaNoWet;
74c4762a1bSJed Brown Field *uold;
75c4762a1bSJed Brown PetscScalar *Kappa;
76c4762a1bSJed Brown PetscInt i;
77c4762a1bSJed Brown
78c4762a1bSJed Brown PetscFunctionBeginUser;
799566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(user->cda, &L));
80c4762a1bSJed Brown
819566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(info->da, user->uold, &uold));
829566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(user->cda, user->Kappa, &Kappa));
83c4762a1bSJed Brown /* Compute residual over the locally owned part of the grid */
84c4762a1bSJed Brown for (i = info->xs; i < info->xs + info->xm; ++i) {
85c4762a1bSJed Brown if (i == 0) {
86c4762a1bSJed Brown f[i].s = u[i].s - user->sl;
87c4762a1bSJed Brown f[i].v = u[i].v - user->vl;
88c4762a1bSJed Brown f[i].p = u[i].p - user->pl;
89c4762a1bSJed Brown } else {
90c4762a1bSJed Brown PetscScalar K = 2 * dx / (dx / Kappa[i] + dx / Kappa[i - 1]);
91c4762a1bSJed Brown PetscReal lambdaWet = kappaWet * PetscRealPart(PetscPowScalar(u[i].s, alpha));
92c4762a1bSJed Brown PetscReal lambda = lambdaWet + kappaNoWet * PetscRealPart(PetscPowScalar(1. - u[i].s, beta));
93c4762a1bSJed Brown PetscReal lambdaWetL = kappaWet * PetscRealPart(PetscPowScalar(u[i - 1].s, alpha));
94c4762a1bSJed Brown PetscReal lambdaL = lambdaWetL + kappaNoWet * PetscRealPart(PetscPowScalar(1. - u[i - 1].s, beta));
95c4762a1bSJed Brown
96c4762a1bSJed Brown f[i].s = phi * (u[i].s - uold[i].s) + (dt / dx) * ((lambdaWet / lambda) * u[i].v - (lambdaWetL / lambdaL) * u[i - 1].v);
97c4762a1bSJed Brown
98c4762a1bSJed Brown f[i].v = u[i].v + K * lambda * (u[i].p - u[i - 1].p) / dx;
99c4762a1bSJed Brown
100c4762a1bSJed Brown /*pxx = (2.0*u[i].p - u[i-1].p - u[i+1].p)/dx;*/
101c4762a1bSJed Brown f[i].p = u[i].v - u[i - 1].v;
102c4762a1bSJed Brown }
103c4762a1bSJed Brown }
1049566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(info->da, user->uold, &uold));
1059566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(user->cda, user->Kappa, &Kappa));
1069566063dSJacob Faibussowitsch /* PetscCall(PetscLogFlops(11.0*info->ym*info->xm)); */
107c4762a1bSJed Brown
1089566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(user->cda, &L));
109*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
110c4762a1bSJed Brown }
111c4762a1bSJed Brown
main(int argc,char ** argv)112d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
113d71ae5a4SJacob Faibussowitsch {
114c4762a1bSJed Brown SNES snes; /* nonlinear solver */
115c4762a1bSJed Brown DM da; /* grid */
116c4762a1bSJed Brown Vec u; /* solution vector */
117c4762a1bSJed Brown AppCtx user; /* user-defined work context */
118c4762a1bSJed Brown PetscReal t = 0.0; /* time */
119c4762a1bSJed Brown PetscInt n;
120c4762a1bSJed Brown
121327415f7SBarry Smith PetscFunctionBeginUser;
1229566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help));
123c4762a1bSJed Brown /* Create solver */
1249566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
125c4762a1bSJed Brown /* Create mesh */
1269566063dSJacob Faibussowitsch PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, 4, 3, 1, NULL, &da));
1279566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da));
1289566063dSJacob Faibussowitsch PetscCall(DMSetUp(da));
1299566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(da, &user));
1309566063dSJacob Faibussowitsch PetscCall(SNESSetDM(snes, da));
131c4762a1bSJed Brown /* Create coefficient */
1329566063dSJacob Faibussowitsch PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, 4, 1, 1, NULL, &user.cda));
1339566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(user.cda));
1349566063dSJacob Faibussowitsch PetscCall(DMSetUp(user.cda));
1359566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(user.cda, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0));
1369566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(user.cda, &user.Kappa));
1379566063dSJacob Faibussowitsch PetscCall(FormPermeability(user.cda, user.Kappa, &user));
138c4762a1bSJed Brown /* Setup Problem */
1399566063dSJacob Faibussowitsch PetscCall(DMDASNESSetFunctionLocal(da, INSERT_VALUES, (PetscErrorCode (*)(DMDALocalInfo *, void *, void *, void *))FormFunctionLocal, &user));
1409566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(da, &u));
1419566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(da, &user.uold));
142c4762a1bSJed Brown
143c4762a1bSJed Brown user.sl = 1.0;
144c4762a1bSJed Brown user.vl = 0.1;
145c4762a1bSJed Brown user.pl = 1.0;
146c4762a1bSJed Brown user.phi = 1.0;
147c4762a1bSJed Brown
148c4762a1bSJed Brown user.kappaWet = 1.0;
149c4762a1bSJed Brown user.kappaNoWet = 0.3;
150c4762a1bSJed Brown
151c4762a1bSJed Brown /* Time Loop */
152c4762a1bSJed Brown user.dt = 0.1;
153c4762a1bSJed Brown for (n = 0; n < 100; ++n, t += user.dt) {
1549566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Starting time %g\n", (double)t));
1559566063dSJacob Faibussowitsch PetscCall(VecView(u, PETSC_VIEWER_DRAW_WORLD));
156c4762a1bSJed Brown /* Solve */
1579566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes));
1589566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, NULL, u));
159c4762a1bSJed Brown /* Update */
1609566063dSJacob Faibussowitsch PetscCall(VecCopy(u, user.uold));
161c4762a1bSJed Brown
1629566063dSJacob Faibussowitsch PetscCall(VecView(u, PETSC_VIEWER_DRAW_WORLD));
163c4762a1bSJed Brown }
164c4762a1bSJed Brown /* Cleanup */
1659566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(da, &u));
1669566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(da, &user.uold));
1679566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(user.cda, &user.Kappa));
1689566063dSJacob Faibussowitsch PetscCall(DMDestroy(&user.cda));
1699566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da));
1709566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes));
1719566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
172b122ec5aSJacob Faibussowitsch return 0;
173c4762a1bSJed Brown }
174c4762a1bSJed Brown
175c4762a1bSJed Brown /*TEST
176c4762a1bSJed Brown
177c4762a1bSJed Brown test:
178c4762a1bSJed Brown suffix: 0
179c4762a1bSJed Brown requires: !single
180c4762a1bSJed Brown args: -snes_converged_reason -snes_monitor_short
181c4762a1bSJed Brown
182c4762a1bSJed Brown TEST*/
183