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 */ 33c4762a1bSJed Brown PetscErrorCode FormPermeability(DM da, Vec Kappa, AppCtx *user) 34c4762a1bSJed Brown { 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; 42*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinateDM(da, &cda)); 43*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinates(da, &c)); 44*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetCorners(da, &xs,NULL,NULL, &xm,NULL,NULL)); 45*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArray(da, Kappa, &K)); 46*5f80ce2aSJacob Faibussowitsch CHKERRQ(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 } 56*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(da, Kappa, &K)); 57*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(cda, c, &coords)); 58c4762a1bSJed Brown PetscFunctionReturn(0); 59c4762a1bSJed Brown } 60c4762a1bSJed Brown 61c4762a1bSJed Brown /* 62c4762a1bSJed Brown FormFunctionLocal - Evaluates nonlinear residual, F(x) on local process patch 63c4762a1bSJed Brown */ 64c4762a1bSJed Brown PetscErrorCode FormFunctionLocal(DMDALocalInfo *info, Field *u, Field *f, AppCtx *user) 65c4762a1bSJed Brown { 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; 79*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetGlobalVector(user->cda, &L)); 80c4762a1bSJed Brown 81*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArray(info->da, user->uold, &uold)); 82*5f80ce2aSJacob Faibussowitsch CHKERRQ(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 } 104*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(info->da, user->uold, &uold)); 105*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(user->cda, user->Kappa, &Kappa)); 106*5f80ce2aSJacob Faibussowitsch /* CHKERRQ(PetscLogFlops(11.0*info->ym*info->xm)); */ 107c4762a1bSJed Brown 108*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreGlobalVector(user->cda, &L)); 109c4762a1bSJed Brown PetscFunctionReturn(0); 110c4762a1bSJed Brown } 111c4762a1bSJed Brown 112c4762a1bSJed Brown int main(int argc, char **argv) 113c4762a1bSJed Brown { 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 PetscErrorCode ierr; 120c4762a1bSJed Brown PetscInt n; 121c4762a1bSJed Brown 122c4762a1bSJed Brown ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr; 123c4762a1bSJed Brown /* Create solver */ 124*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESCreate(PETSC_COMM_WORLD, &snes)); 125c4762a1bSJed Brown /* Create mesh */ 126*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,4,3,1,NULL,&da)); 127*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(da)); 128*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(da)); 129*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetApplicationContext(da, &user)); 130*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetDM(snes, da)); 131c4762a1bSJed Brown /* Create coefficient */ 132*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,4,1,1,NULL,&user.cda)); 133*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(user.cda)); 134*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(user.cda)); 135*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDASetUniformCoordinates(user.cda, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0)); 136*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetGlobalVector(user.cda, &user.Kappa)); 137*5f80ce2aSJacob Faibussowitsch CHKERRQ(FormPermeability(user.cda, user.Kappa, &user)); 138c4762a1bSJed Brown /* Setup Problem */ 139*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDASNESSetFunctionLocal(da,INSERT_VALUES,(PetscErrorCode (*)(DMDALocalInfo*,void*,void*,void*))FormFunctionLocal,&user)); 140*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetGlobalVector(da, &u)); 141*5f80ce2aSJacob Faibussowitsch CHKERRQ(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) { 154*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "Starting time %g\n", (double)t)); 155*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(u, PETSC_VIEWER_DRAW_WORLD)); 156c4762a1bSJed Brown /* Solve */ 157*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetFromOptions(snes)); 158*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSolve(snes, NULL, u)); 159c4762a1bSJed Brown /* Update */ 160*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(u, user.uold)); 161c4762a1bSJed Brown 162*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(u, PETSC_VIEWER_DRAW_WORLD)); 163c4762a1bSJed Brown } 164c4762a1bSJed Brown /* Cleanup */ 165*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreGlobalVector(da, &u)); 166*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreGlobalVector(da, &user.uold)); 167*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreGlobalVector(user.cda, &user.Kappa)); 168*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&user.cda)); 169*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&da)); 170*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESDestroy(&snes)); 171c4762a1bSJed Brown ierr = PetscFinalize(); 172c4762a1bSJed Brown return ierr; 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