1c4762a1bSJed Brown static char help[] = "Bratu nonlinear PDE in 2d.\n\ 2c4762a1bSJed Brown We solve the Bratu (SFI - solid fuel ignition) problem in a 2D rectangular\n\ 3c4762a1bSJed Brown domain, using distributed arrays (DMDAs) to partition the parallel grid.\n\ 4c4762a1bSJed Brown The command line options include:\n\ 5c4762a1bSJed Brown -par <parameter>, where <parameter> indicates the problem's nonlinearity\n\ 6c4762a1bSJed Brown problem SFI: <parameter> = Bratu parameter (0 <= par <= 6.81)\n\n\ 7c4762a1bSJed Brown -m_par/n_par <parameter>, where <parameter> indicates an integer\n \ 8c4762a1bSJed Brown that MMS3 will be evaluated with 2^m_par, 2^n_par"; 9c4762a1bSJed Brown 10c4762a1bSJed Brown /*T 11c4762a1bSJed Brown Concepts: SNES^parallel Bratu example 12c4762a1bSJed Brown Concepts: DMDA^using distributed arrays; 13c4762a1bSJed Brown Concepts: IS coloirng types; 14c4762a1bSJed Brown Processors: n 15c4762a1bSJed Brown T*/ 16c4762a1bSJed Brown 17c4762a1bSJed Brown 18c4762a1bSJed Brown 19c4762a1bSJed Brown /* ------------------------------------------------------------------------ 20c4762a1bSJed Brown 21c4762a1bSJed Brown Solid Fuel Ignition (SFI) problem. This problem is modeled by 22c4762a1bSJed Brown the partial differential equation 23c4762a1bSJed Brown 24c4762a1bSJed Brown -Laplacian u - lambda*exp(u) = 0, 0 < x,y < 1, 25c4762a1bSJed Brown 26c4762a1bSJed Brown with boundary conditions 27c4762a1bSJed Brown 28c4762a1bSJed Brown u = 0 for x = 0, x = 1, y = 0, y = 1. 29c4762a1bSJed Brown 30c4762a1bSJed Brown A finite difference approximation with the usual 5-point stencil 31c4762a1bSJed Brown is used to discretize the boundary value problem to obtain a nonlinear 32c4762a1bSJed Brown system of equations. 33c4762a1bSJed Brown 34c4762a1bSJed Brown 35c4762a1bSJed Brown This example shows how geometric multigrid can be run transparently with a nonlinear solver so long 36c4762a1bSJed Brown as SNESSetDM() is provided. Example usage 37c4762a1bSJed Brown 38c4762a1bSJed Brown ./ex5 -pc_type mg -ksp_monitor -snes_view -pc_mg_levels 3 -pc_mg_galerkin pmat -da_grid_x 17 -da_grid_y 17 39c4762a1bSJed Brown -mg_levels_ksp_monitor -snes_monitor -mg_levels_pc_type sor -pc_mg_type full 40c4762a1bSJed Brown 41c4762a1bSJed Brown or to run with grid sequencing on the nonlinear problem (note that you do not need to provide the number of 42c4762a1bSJed Brown multigrid levels, it will be determined automatically based on the number of refinements done) 43c4762a1bSJed Brown 44c4762a1bSJed Brown ./ex5 -pc_type mg -ksp_monitor -snes_view -pc_mg_galerkin pmat -snes_grid_sequence 3 45c4762a1bSJed Brown -mg_levels_ksp_monitor -snes_monitor -mg_levels_pc_type sor -pc_mg_type full 46c4762a1bSJed Brown 47c4762a1bSJed Brown ------------------------------------------------------------------------- */ 48c4762a1bSJed Brown 49c4762a1bSJed Brown /* 50c4762a1bSJed Brown Include "petscdmda.h" so that we can use distributed arrays (DMDAs). 51c4762a1bSJed Brown Include "petscsnes.h" so that we can use SNES solvers. Note that this 52c4762a1bSJed Brown */ 53c4762a1bSJed Brown #include <petscdm.h> 54c4762a1bSJed Brown #include <petscdmda.h> 55c4762a1bSJed Brown #include <petscsnes.h> 56c4762a1bSJed Brown #include <petscmatlab.h> 57c4762a1bSJed Brown #include <petsc/private/snesimpl.h> /* For SNES_Solve event */ 58c4762a1bSJed Brown 59c4762a1bSJed Brown /* 60c4762a1bSJed Brown User-defined application context - contains data needed by the 61c4762a1bSJed Brown application-provided call-back routines, FormJacobianLocal() and 62c4762a1bSJed Brown FormFunctionLocal(). 63c4762a1bSJed Brown */ 64c4762a1bSJed Brown typedef struct AppCtx AppCtx; 65c4762a1bSJed Brown struct AppCtx { 66c4762a1bSJed Brown PetscReal param; /* test problem parameter */ 67c4762a1bSJed Brown PetscInt m,n; /* MMS3 parameters */ 68c4762a1bSJed Brown PetscErrorCode (*mms_solution)(AppCtx*,const DMDACoor2d*,PetscScalar*); 69c4762a1bSJed Brown PetscErrorCode (*mms_forcing)(AppCtx*,const DMDACoor2d*,PetscScalar*); 70c4762a1bSJed Brown }; 71c4762a1bSJed Brown 72c4762a1bSJed Brown /* 73c4762a1bSJed Brown User-defined routines 74c4762a1bSJed Brown */ 75c4762a1bSJed Brown extern PetscErrorCode FormInitialGuess(DM,AppCtx*,Vec); 76c4762a1bSJed Brown extern PetscErrorCode FormFunctionLocal(DMDALocalInfo*,PetscScalar**,PetscScalar**,AppCtx*); 77c4762a1bSJed Brown extern PetscErrorCode FormExactSolution(DM,AppCtx*,Vec); 78c4762a1bSJed Brown extern PetscErrorCode ZeroBCSolution(AppCtx*,const DMDACoor2d*,PetscScalar*); 79c4762a1bSJed Brown extern PetscErrorCode MMSSolution1(AppCtx*,const DMDACoor2d*,PetscScalar*); 80c4762a1bSJed Brown extern PetscErrorCode MMSForcing1(AppCtx*,const DMDACoor2d*,PetscScalar*); 81c4762a1bSJed Brown extern PetscErrorCode MMSSolution2(AppCtx*,const DMDACoor2d*,PetscScalar*); 82c4762a1bSJed Brown extern PetscErrorCode MMSForcing2(AppCtx*,const DMDACoor2d*,PetscScalar*); 83c4762a1bSJed Brown extern PetscErrorCode MMSSolution3(AppCtx*,const DMDACoor2d*,PetscScalar*); 84c4762a1bSJed Brown extern PetscErrorCode MMSForcing3(AppCtx*,const DMDACoor2d*,PetscScalar*); 85c4762a1bSJed Brown extern PetscErrorCode MMSSolution4(AppCtx*,const DMDACoor2d*,PetscScalar*); 86c4762a1bSJed Brown extern PetscErrorCode MMSForcing4(AppCtx*,const DMDACoor2d*,PetscScalar*); 87c4762a1bSJed Brown extern PetscErrorCode FormJacobianLocal(DMDALocalInfo*,PetscScalar**,Mat,Mat,AppCtx*); 88c4762a1bSJed Brown extern PetscErrorCode FormObjectiveLocal(DMDALocalInfo*,PetscScalar**,PetscReal*,AppCtx*); 89c4762a1bSJed Brown extern PetscErrorCode FormFunctionMatlab(SNES,Vec,Vec,void*); 90c4762a1bSJed Brown extern PetscErrorCode NonlinearGS(SNES,Vec,Vec,void*); 91c4762a1bSJed Brown 92c4762a1bSJed Brown int main(int argc,char **argv) 93c4762a1bSJed Brown { 94c4762a1bSJed Brown SNES snes; /* nonlinear solver */ 95c4762a1bSJed Brown Vec x; /* solution vector */ 96c4762a1bSJed Brown AppCtx user; /* user-defined work context */ 97c4762a1bSJed Brown PetscInt its; /* iterations for convergence */ 98c4762a1bSJed Brown PetscErrorCode ierr; 99c4762a1bSJed Brown PetscReal bratu_lambda_max = 6.81; 100c4762a1bSJed Brown PetscReal bratu_lambda_min = 0.; 101c4762a1bSJed Brown PetscInt MMS = 0; 102c4762a1bSJed Brown PetscBool flg = PETSC_FALSE; 103c4762a1bSJed Brown DM da; 104c4762a1bSJed Brown Vec r = NULL; 105c4762a1bSJed Brown KSP ksp; 106c4762a1bSJed Brown PetscInt lits,slits; 107c4762a1bSJed Brown 108c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 109c4762a1bSJed Brown Initialize program 110c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 111c4762a1bSJed Brown 112c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 113c4762a1bSJed Brown 114c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 115c4762a1bSJed Brown Initialize problem parameters 116c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 117c4762a1bSJed Brown user.param = 6.0; 118c4762a1bSJed Brown ierr = PetscOptionsGetReal(NULL,NULL,"-par",&user.param,NULL);CHKERRQ(ierr); 119c4762a1bSJed Brown if (user.param > bratu_lambda_max || user.param < bratu_lambda_min) SETERRQ3(PETSC_COMM_SELF,1,"Lambda, %g, is out of range, [%g, %g]", user.param, bratu_lambda_min, bratu_lambda_max); 120c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-mms",&MMS,NULL);CHKERRQ(ierr); 121c4762a1bSJed Brown if (MMS == 3) { 122c4762a1bSJed Brown PetscInt mPar = 2, nPar = 1; 123c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-m_par",&mPar,NULL);CHKERRQ(ierr); 124c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-n_par",&nPar,NULL);CHKERRQ(ierr); 125c4762a1bSJed Brown user.m = PetscPowInt(2,mPar); 126c4762a1bSJed Brown user.n = PetscPowInt(2,nPar); 127c4762a1bSJed Brown } 128c4762a1bSJed Brown 129c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 130c4762a1bSJed Brown Create nonlinear solver context 131c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 132c4762a1bSJed Brown ierr = SNESCreate(PETSC_COMM_WORLD,&snes);CHKERRQ(ierr); 133c4762a1bSJed Brown ierr = SNESSetCountersReset(snes,PETSC_FALSE);CHKERRQ(ierr); 134c4762a1bSJed Brown ierr = SNESSetNGS(snes, NonlinearGS, NULL);CHKERRQ(ierr); 135c4762a1bSJed Brown 136c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 137c4762a1bSJed Brown Create distributed array (DMDA) to manage parallel grid and vectors 138c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 139c4762a1bSJed Brown ierr = DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,4,4,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&da);CHKERRQ(ierr); 140c4762a1bSJed Brown ierr = DMSetFromOptions(da);CHKERRQ(ierr); 141c4762a1bSJed Brown ierr = DMSetUp(da);CHKERRQ(ierr); 142c4762a1bSJed Brown ierr = DMDASetUniformCoordinates(da, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0);CHKERRQ(ierr); 143c4762a1bSJed Brown ierr = DMSetApplicationContext(da,&user);CHKERRQ(ierr); 144c4762a1bSJed Brown ierr = SNESSetDM(snes,da);CHKERRQ(ierr); 145c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 146c4762a1bSJed Brown Extract global vectors from DMDA; then duplicate for remaining 147c4762a1bSJed Brown vectors that are the same types 148c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 149c4762a1bSJed Brown ierr = DMCreateGlobalVector(da,&x);CHKERRQ(ierr); 150c4762a1bSJed Brown 151c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 152c4762a1bSJed Brown Set local function evaluation routine 153c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 154c4762a1bSJed Brown user.mms_solution = ZeroBCSolution; 155c4762a1bSJed Brown switch (MMS) { 156c4762a1bSJed Brown case 0: user.mms_solution = NULL; user.mms_forcing = NULL; CHKERRQ(ierr); 157c4762a1bSJed Brown case 1: user.mms_solution = MMSSolution1; user.mms_forcing = MMSForcing1; break; 158c4762a1bSJed Brown case 2: user.mms_solution = MMSSolution2; user.mms_forcing = MMSForcing2; break; 159c4762a1bSJed Brown case 3: user.mms_solution = MMSSolution3; user.mms_forcing = MMSForcing3; break; 160c4762a1bSJed Brown case 4: user.mms_solution = MMSSolution4; user.mms_forcing = MMSForcing4; break; 161c4762a1bSJed Brown default: SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_USER,"Unknown MMS type %d",MMS); 162c4762a1bSJed Brown } 163c4762a1bSJed Brown ierr = DMDASNESSetFunctionLocal(da,INSERT_VALUES,(DMDASNESFunction)FormFunctionLocal,&user);CHKERRQ(ierr); 164c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-fd",&flg,NULL);CHKERRQ(ierr); 165c4762a1bSJed Brown if (!flg) { 166c4762a1bSJed Brown ierr = DMDASNESSetJacobianLocal(da,(DMDASNESJacobian)FormJacobianLocal,&user);CHKERRQ(ierr); 167c4762a1bSJed Brown } 168c4762a1bSJed Brown 169c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-obj",&flg,NULL);CHKERRQ(ierr); 170c4762a1bSJed Brown if (flg) { 171c4762a1bSJed Brown ierr = DMDASNESSetObjectiveLocal(da,(DMDASNESObjective)FormObjectiveLocal,&user);CHKERRQ(ierr); 172c4762a1bSJed Brown } 173c4762a1bSJed Brown 174*4e1ad211SJed Brown if (PetscDefined(HAVE_MATLAB_ENGINE)) { 175*4e1ad211SJed Brown PetscBool matlab_function = PETSC_FALSE; 176c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-matlab_function",&matlab_function,0);CHKERRQ(ierr); 177c4762a1bSJed Brown if (matlab_function) { 178c4762a1bSJed Brown ierr = VecDuplicate(x,&r);CHKERRQ(ierr); 179c4762a1bSJed Brown ierr = SNESSetFunction(snes,r,FormFunctionMatlab,&user);CHKERRQ(ierr); 180c4762a1bSJed Brown } 181*4e1ad211SJed Brown } 182c4762a1bSJed Brown 183c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 184c4762a1bSJed Brown Customize nonlinear solver; set runtime options 185c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 186c4762a1bSJed Brown ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 187c4762a1bSJed Brown 188c4762a1bSJed Brown ierr = FormInitialGuess(da,&user,x);CHKERRQ(ierr); 189c4762a1bSJed Brown 190c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 191c4762a1bSJed Brown Solve nonlinear system 192c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 193c4762a1bSJed Brown ierr = SNESSolve(snes,NULL,x);CHKERRQ(ierr); 194c4762a1bSJed Brown ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr); 195c4762a1bSJed Brown 196c4762a1bSJed Brown ierr = SNESGetLinearSolveIterations(snes,&slits);CHKERRQ(ierr); 197c4762a1bSJed Brown ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 198c4762a1bSJed Brown ierr = KSPGetTotalIterations(ksp,&lits);CHKERRQ(ierr); 199c4762a1bSJed Brown if (lits != slits) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Number of total linear iterations reported by SNES %D does not match reported by KSP %D",slits,lits); 200c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 201c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 202c4762a1bSJed Brown 203c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 204c4762a1bSJed Brown If using MMS, check the l_2 error 205c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 206c4762a1bSJed Brown if (MMS) { 207c4762a1bSJed Brown Vec e; 208c4762a1bSJed Brown PetscReal errorl2, errorinf; 209c4762a1bSJed Brown PetscInt N; 210c4762a1bSJed Brown 211c4762a1bSJed Brown ierr = VecDuplicate(x, &e);CHKERRQ(ierr); 212c4762a1bSJed Brown ierr = PetscObjectViewFromOptions((PetscObject) x, NULL, "-sol_view");CHKERRQ(ierr); 213c4762a1bSJed Brown ierr = FormExactSolution(da, &user, e);CHKERRQ(ierr); 214c4762a1bSJed Brown ierr = PetscObjectViewFromOptions((PetscObject) e, NULL, "-exact_view");CHKERRQ(ierr); 215c4762a1bSJed Brown ierr = VecAXPY(e, -1.0, x);CHKERRQ(ierr); 216c4762a1bSJed Brown ierr = PetscObjectViewFromOptions((PetscObject) e, NULL, "-error_view");CHKERRQ(ierr); 217c4762a1bSJed Brown ierr = VecNorm(e, NORM_2, &errorl2);CHKERRQ(ierr); 218c4762a1bSJed Brown ierr = VecNorm(e, NORM_INFINITY, &errorinf);CHKERRQ(ierr); 219c4762a1bSJed Brown ierr = VecGetSize(e, &N);CHKERRQ(ierr); 220c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD, "N: %D error L2 %g inf %g\n", N, (double) errorl2/PetscSqrtReal(N), (double) errorinf);CHKERRQ(ierr); 221c4762a1bSJed Brown ierr = VecDestroy(&e);CHKERRQ(ierr); 222c4762a1bSJed Brown ierr = PetscLogEventSetDof(SNES_Solve, 0, N);CHKERRQ(ierr); 223c4762a1bSJed Brown ierr = PetscLogEventSetError(SNES_Solve, 0, errorl2/PetscSqrtReal(N));CHKERRQ(ierr); 224c4762a1bSJed Brown } 225c4762a1bSJed Brown 226c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 227c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 228c4762a1bSJed Brown are no longer needed. 229c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 230c4762a1bSJed Brown ierr = VecDestroy(&r);CHKERRQ(ierr); 231c4762a1bSJed Brown ierr = VecDestroy(&x);CHKERRQ(ierr); 232c4762a1bSJed Brown ierr = SNESDestroy(&snes);CHKERRQ(ierr); 233c4762a1bSJed Brown ierr = DMDestroy(&da);CHKERRQ(ierr); 234c4762a1bSJed Brown ierr = PetscFinalize(); 235c4762a1bSJed Brown return ierr; 236c4762a1bSJed Brown } 237c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 238c4762a1bSJed Brown /* 239c4762a1bSJed Brown FormInitialGuess - Forms initial approximation. 240c4762a1bSJed Brown 241c4762a1bSJed Brown Input Parameters: 242c4762a1bSJed Brown da - The DM 243c4762a1bSJed Brown user - user-defined application context 244c4762a1bSJed Brown 245c4762a1bSJed Brown Output Parameter: 246c4762a1bSJed Brown X - vector 247c4762a1bSJed Brown */ 248c4762a1bSJed Brown PetscErrorCode FormInitialGuess(DM da,AppCtx *user,Vec X) 249c4762a1bSJed Brown { 250c4762a1bSJed Brown PetscInt i,j,Mx,My,xs,ys,xm,ym; 251c4762a1bSJed Brown PetscErrorCode ierr; 252c4762a1bSJed Brown PetscReal lambda,temp1,temp,hx,hy; 253c4762a1bSJed Brown PetscScalar **x; 254c4762a1bSJed Brown 255c4762a1bSJed Brown PetscFunctionBeginUser; 256c4762a1bSJed Brown ierr = DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);CHKERRQ(ierr); 257c4762a1bSJed Brown 258c4762a1bSJed Brown lambda = user->param; 259c4762a1bSJed Brown hx = 1.0/(PetscReal)(Mx-1); 260c4762a1bSJed Brown hy = 1.0/(PetscReal)(My-1); 261c4762a1bSJed Brown temp1 = lambda/(lambda + 1.0); 262c4762a1bSJed Brown 263c4762a1bSJed Brown /* 264c4762a1bSJed Brown Get a pointer to vector data. 265c4762a1bSJed Brown - For default PETSc vectors, VecGetArray() returns a pointer to 266c4762a1bSJed Brown the data array. Otherwise, the routine is implementation dependent. 267c4762a1bSJed Brown - You MUST call VecRestoreArray() when you no longer need access to 268c4762a1bSJed Brown the array. 269c4762a1bSJed Brown */ 270c4762a1bSJed Brown ierr = DMDAVecGetArray(da,X,&x);CHKERRQ(ierr); 271c4762a1bSJed Brown 272c4762a1bSJed Brown /* 273c4762a1bSJed Brown Get local grid boundaries (for 2-dimensional DMDA): 274c4762a1bSJed Brown xs, ys - starting grid indices (no ghost points) 275c4762a1bSJed Brown xm, ym - widths of local grid (no ghost points) 276c4762a1bSJed Brown 277c4762a1bSJed Brown */ 278c4762a1bSJed Brown ierr = DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr); 279c4762a1bSJed Brown 280c4762a1bSJed Brown /* 281c4762a1bSJed Brown Compute initial guess over the locally owned part of the grid 282c4762a1bSJed Brown */ 283c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 284c4762a1bSJed Brown temp = (PetscReal)(PetscMin(j,My-j-1))*hy; 285c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 286c4762a1bSJed Brown if (i == 0 || j == 0 || i == Mx-1 || j == My-1) { 287c4762a1bSJed Brown /* boundary conditions are all zero Dirichlet */ 288c4762a1bSJed Brown x[j][i] = 0.0; 289c4762a1bSJed Brown } else { 290c4762a1bSJed Brown x[j][i] = temp1*PetscSqrtReal(PetscMin((PetscReal)(PetscMin(i,Mx-i-1))*hx,temp)); 291c4762a1bSJed Brown } 292c4762a1bSJed Brown } 293c4762a1bSJed Brown } 294c4762a1bSJed Brown 295c4762a1bSJed Brown /* 296c4762a1bSJed Brown Restore vector 297c4762a1bSJed Brown */ 298c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da,X,&x);CHKERRQ(ierr); 299c4762a1bSJed Brown PetscFunctionReturn(0); 300c4762a1bSJed Brown } 301c4762a1bSJed Brown 302c4762a1bSJed Brown /* 303c4762a1bSJed Brown FormExactSolution - Forms MMS solution 304c4762a1bSJed Brown 305c4762a1bSJed Brown Input Parameters: 306c4762a1bSJed Brown da - The DM 307c4762a1bSJed Brown user - user-defined application context 308c4762a1bSJed Brown 309c4762a1bSJed Brown Output Parameter: 310c4762a1bSJed Brown X - vector 311c4762a1bSJed Brown */ 312c4762a1bSJed Brown PetscErrorCode FormExactSolution(DM da, AppCtx *user, Vec U) 313c4762a1bSJed Brown { 314c4762a1bSJed Brown DM coordDA; 315c4762a1bSJed Brown Vec coordinates; 316c4762a1bSJed Brown DMDACoor2d **coords; 317c4762a1bSJed Brown PetscScalar **u; 318c4762a1bSJed Brown PetscInt xs, ys, xm, ym, i, j; 319c4762a1bSJed Brown PetscErrorCode ierr; 320c4762a1bSJed Brown 321c4762a1bSJed Brown PetscFunctionBeginUser; 322c4762a1bSJed Brown ierr = DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL);CHKERRQ(ierr); 323c4762a1bSJed Brown ierr = DMGetCoordinateDM(da, &coordDA);CHKERRQ(ierr); 324c4762a1bSJed Brown ierr = DMGetCoordinates(da, &coordinates);CHKERRQ(ierr); 325c4762a1bSJed Brown ierr = DMDAVecGetArray(coordDA, coordinates, &coords);CHKERRQ(ierr); 326c4762a1bSJed Brown ierr = DMDAVecGetArray(da, U, &u);CHKERRQ(ierr); 327c4762a1bSJed Brown for (j = ys; j < ys+ym; ++j) { 328c4762a1bSJed Brown for (i = xs; i < xs+xm; ++i) { 329c4762a1bSJed Brown user->mms_solution(user,&coords[j][i],&u[j][i]); 330c4762a1bSJed Brown } 331c4762a1bSJed Brown } 332c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da, U, &u);CHKERRQ(ierr); 333c4762a1bSJed Brown ierr = DMDAVecRestoreArray(coordDA, coordinates, &coords);CHKERRQ(ierr); 334c4762a1bSJed Brown PetscFunctionReturn(0); 335c4762a1bSJed Brown } 336c4762a1bSJed Brown 337c4762a1bSJed Brown PetscErrorCode ZeroBCSolution(AppCtx *user,const DMDACoor2d *c,PetscScalar *u) 338c4762a1bSJed Brown { 339c4762a1bSJed Brown u[0] = 0.; 340c4762a1bSJed Brown return 0; 341c4762a1bSJed Brown } 342c4762a1bSJed Brown 343c4762a1bSJed Brown /* The functions below evaluate the MMS solution u(x,y) and associated forcing 344c4762a1bSJed Brown 345c4762a1bSJed Brown f(x,y) = -u_xx - u_yy - lambda exp(u) 346c4762a1bSJed Brown 347c4762a1bSJed Brown such that u(x,y) is an exact solution with f(x,y) as the right hand side forcing term. 348c4762a1bSJed Brown */ 349c4762a1bSJed Brown PetscErrorCode MMSSolution1(AppCtx *user,const DMDACoor2d *c,PetscScalar *u) 350c4762a1bSJed Brown { 351c4762a1bSJed Brown PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y); 352c4762a1bSJed Brown u[0] = x*(1 - x)*y*(1 - y); 353c4762a1bSJed Brown PetscLogFlops(5); 354c4762a1bSJed Brown return 0; 355c4762a1bSJed Brown } 356c4762a1bSJed Brown PetscErrorCode MMSForcing1(AppCtx *user,const DMDACoor2d *c,PetscScalar *f) 357c4762a1bSJed Brown { 358c4762a1bSJed Brown PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y); 359c4762a1bSJed Brown f[0] = 2*x*(1 - x) + 2*y*(1 - y) - user->param*PetscExpReal(x*(1 - x)*y*(1 - y)); 360c4762a1bSJed Brown return 0; 361c4762a1bSJed Brown } 362c4762a1bSJed Brown 363c4762a1bSJed Brown PetscErrorCode MMSSolution2(AppCtx *user,const DMDACoor2d *c,PetscScalar *u) 364c4762a1bSJed Brown { 365c4762a1bSJed Brown PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y); 366c4762a1bSJed Brown u[0] = PetscSinReal(PETSC_PI*x)*PetscSinReal(PETSC_PI*y); 367c4762a1bSJed Brown PetscLogFlops(5); 368c4762a1bSJed Brown return 0; 369c4762a1bSJed Brown } 370c4762a1bSJed Brown PetscErrorCode MMSForcing2(AppCtx *user,const DMDACoor2d *c,PetscScalar *f) 371c4762a1bSJed Brown { 372c4762a1bSJed Brown PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y); 373c4762a1bSJed Brown f[0] = 2*PetscSqr(PETSC_PI)*PetscSinReal(PETSC_PI*x)*PetscSinReal(PETSC_PI*y) - user->param*PetscExpReal(PetscSinReal(PETSC_PI*x)*PetscSinReal(PETSC_PI*y)); 374c4762a1bSJed Brown return 0; 375c4762a1bSJed Brown } 376c4762a1bSJed Brown 377c4762a1bSJed Brown PetscErrorCode MMSSolution3(AppCtx *user,const DMDACoor2d *c,PetscScalar *u) 378c4762a1bSJed Brown { 379c4762a1bSJed Brown PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y); 380c4762a1bSJed Brown u[0] = PetscSinReal(user->m*PETSC_PI*x*(1-y))*PetscSinReal(user->n*PETSC_PI*y*(1-x)); 381c4762a1bSJed Brown PetscLogFlops(5); 382c4762a1bSJed Brown return 0; 383c4762a1bSJed Brown } 384c4762a1bSJed Brown PetscErrorCode MMSForcing3(AppCtx *user,const DMDACoor2d *c,PetscScalar *f) 385c4762a1bSJed Brown { 386c4762a1bSJed Brown PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y); 387c4762a1bSJed Brown PetscReal m = user->m, n = user->n, lambda = user->param; 388c4762a1bSJed Brown f[0] = (-(PetscExpReal(PetscSinReal(m*PETSC_PI*x*(1 - y))*PetscSinReal(n*PETSC_PI*(1 - x)*y))*lambda) 389c4762a1bSJed Brown + PetscSqr(PETSC_PI)*(-2*m*n*((-1 + x)*x + (-1 + y)*y)*PetscCosReal(m*PETSC_PI*x*(-1 + y))*PetscCosReal(n*PETSC_PI*(-1 + x)*y) 390c4762a1bSJed Brown + (PetscSqr(m)*(PetscSqr(x) + PetscSqr(-1 + y)) + PetscSqr(n)*(PetscSqr(-1 + x) + PetscSqr(y))) 391c4762a1bSJed Brown *PetscSinReal(m*PETSC_PI*x*(-1 + y))*PetscSinReal(n*PETSC_PI*(-1 + x)*y))); 392c4762a1bSJed Brown return 0; 393c4762a1bSJed Brown } 394c4762a1bSJed Brown 395c4762a1bSJed Brown PetscErrorCode MMSSolution4(AppCtx *user,const DMDACoor2d *c,PetscScalar *u) 396c4762a1bSJed Brown { 397c4762a1bSJed Brown const PetscReal Lx = 1.,Ly = 1.; 398c4762a1bSJed Brown PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y); 399c4762a1bSJed Brown u[0] = (PetscPowReal(x,4)-PetscSqr(Lx)*PetscSqr(x))*(PetscPowReal(y,4)-PetscSqr(Ly)*PetscSqr(y)); 400c4762a1bSJed Brown PetscLogFlops(9); 401c4762a1bSJed Brown return 0; 402c4762a1bSJed Brown } 403c4762a1bSJed Brown PetscErrorCode MMSForcing4(AppCtx *user,const DMDACoor2d *c,PetscScalar *f) 404c4762a1bSJed Brown { 405c4762a1bSJed Brown const PetscReal Lx = 1.,Ly = 1.; 406c4762a1bSJed Brown PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y); 407c4762a1bSJed Brown f[0] = (2*PetscSqr(x)*(PetscSqr(x)-PetscSqr(Lx))*(PetscSqr(Ly)-6*PetscSqr(y)) 408c4762a1bSJed Brown + 2*PetscSqr(y)*(PetscSqr(Lx)-6*PetscSqr(x))*(PetscSqr(y)-PetscSqr(Ly)) 409c4762a1bSJed Brown - user->param*PetscExpReal((PetscPowReal(x,4)-PetscSqr(Lx)*PetscSqr(x))*(PetscPowReal(y,4)-PetscSqr(Ly)*PetscSqr(y)))); 410c4762a1bSJed Brown return 0; 411c4762a1bSJed Brown } 412c4762a1bSJed Brown 413c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 414c4762a1bSJed Brown /* 415c4762a1bSJed Brown FormFunctionLocal - Evaluates nonlinear function, F(x) on local process patch 416c4762a1bSJed Brown 417c4762a1bSJed Brown 418c4762a1bSJed Brown */ 419c4762a1bSJed Brown PetscErrorCode FormFunctionLocal(DMDALocalInfo *info,PetscScalar **x,PetscScalar **f,AppCtx *user) 420c4762a1bSJed Brown { 421c4762a1bSJed Brown PetscErrorCode ierr; 422c4762a1bSJed Brown PetscInt i,j; 423c4762a1bSJed Brown PetscReal lambda,hx,hy,hxdhy,hydhx; 424c4762a1bSJed Brown PetscScalar u,ue,uw,un,us,uxx,uyy,mms_solution,mms_forcing; 425c4762a1bSJed Brown DMDACoor2d c; 426c4762a1bSJed Brown 427c4762a1bSJed Brown PetscFunctionBeginUser; 428c4762a1bSJed Brown lambda = user->param; 429c4762a1bSJed Brown hx = 1.0/(PetscReal)(info->mx-1); 430c4762a1bSJed Brown hy = 1.0/(PetscReal)(info->my-1); 431c4762a1bSJed Brown hxdhy = hx/hy; 432c4762a1bSJed Brown hydhx = hy/hx; 433c4762a1bSJed Brown /* 434c4762a1bSJed Brown Compute function over the locally owned part of the grid 435c4762a1bSJed Brown */ 436c4762a1bSJed Brown for (j=info->ys; j<info->ys+info->ym; j++) { 437c4762a1bSJed Brown for (i=info->xs; i<info->xs+info->xm; i++) { 438c4762a1bSJed Brown if (i == 0 || j == 0 || i == info->mx-1 || j == info->my-1) { 439c4762a1bSJed Brown c.x = i*hx; c.y = j*hy; 440c4762a1bSJed Brown ierr = user->mms_solution(user,&c,&mms_solution);CHKERRQ(ierr); 441c4762a1bSJed Brown f[j][i] = 2.0*(hydhx+hxdhy)*(x[j][i] - mms_solution); 442c4762a1bSJed Brown } else { 443c4762a1bSJed Brown u = x[j][i]; 444c4762a1bSJed Brown uw = x[j][i-1]; 445c4762a1bSJed Brown ue = x[j][i+1]; 446c4762a1bSJed Brown un = x[j-1][i]; 447c4762a1bSJed Brown us = x[j+1][i]; 448c4762a1bSJed Brown 449c4762a1bSJed Brown /* Enforce boundary conditions at neighboring points -- setting these values causes the Jacobian to be symmetric. */ 450c4762a1bSJed Brown if (i-1 == 0) {c.x = (i-1)*hx; c.y = j*hy; ierr = user->mms_solution(user,&c,&uw);CHKERRQ(ierr);} 451c4762a1bSJed Brown if (i+1 == info->mx-1) {c.x = (i+1)*hx; c.y = j*hy; ierr = user->mms_solution(user,&c,&ue);CHKERRQ(ierr);} 452c4762a1bSJed Brown if (j-1 == 0) {c.x = i*hx; c.y = (j-1)*hy; ierr = user->mms_solution(user,&c,&un);CHKERRQ(ierr);} 453c4762a1bSJed Brown if (j+1 == info->my-1) {c.x = i*hx; c.y = (j+1)*hy; ierr = user->mms_solution(user,&c,&us);CHKERRQ(ierr);} 454c4762a1bSJed Brown 455c4762a1bSJed Brown uxx = (2.0*u - uw - ue)*hydhx; 456c4762a1bSJed Brown uyy = (2.0*u - un - us)*hxdhy; 457c4762a1bSJed Brown mms_forcing = 0; 458c4762a1bSJed Brown c.x = i*hx; c.y = j*hy; 459c4762a1bSJed Brown if (user->mms_forcing) {ierr = user->mms_forcing(user,&c,&mms_forcing);CHKERRQ(ierr);} 460c4762a1bSJed Brown f[j][i] = uxx + uyy - hx*hy*(lambda*PetscExpScalar(u) + mms_forcing); 461c4762a1bSJed Brown } 462c4762a1bSJed Brown } 463c4762a1bSJed Brown } 464c4762a1bSJed Brown ierr = PetscLogFlops(11.0*info->ym*info->xm);CHKERRQ(ierr); 465c4762a1bSJed Brown PetscFunctionReturn(0); 466c4762a1bSJed Brown } 467c4762a1bSJed Brown 468c4762a1bSJed Brown /* FormObjectiveLocal - Evaluates nonlinear function, F(x) on local process patch */ 469c4762a1bSJed Brown PetscErrorCode FormObjectiveLocal(DMDALocalInfo *info,PetscScalar **x,PetscReal *obj,AppCtx *user) 470c4762a1bSJed Brown { 471c4762a1bSJed Brown PetscErrorCode ierr; 472c4762a1bSJed Brown PetscInt i,j; 473c4762a1bSJed Brown PetscReal lambda,hx,hy,hxdhy,hydhx,sc,lobj=0; 474c4762a1bSJed Brown PetscScalar u,ue,uw,un,us,uxux,uyuy; 475c4762a1bSJed Brown MPI_Comm comm; 476c4762a1bSJed Brown 477c4762a1bSJed Brown PetscFunctionBeginUser; 478c4762a1bSJed Brown *obj = 0; 479c4762a1bSJed Brown ierr = PetscObjectGetComm((PetscObject)info->da,&comm);CHKERRQ(ierr); 480c4762a1bSJed Brown lambda = user->param; 481c4762a1bSJed Brown hx = 1.0/(PetscReal)(info->mx-1); 482c4762a1bSJed Brown hy = 1.0/(PetscReal)(info->my-1); 483c4762a1bSJed Brown sc = hx*hy*lambda; 484c4762a1bSJed Brown hxdhy = hx/hy; 485c4762a1bSJed Brown hydhx = hy/hx; 486c4762a1bSJed Brown /* 487c4762a1bSJed Brown Compute function over the locally owned part of the grid 488c4762a1bSJed Brown */ 489c4762a1bSJed Brown for (j=info->ys; j<info->ys+info->ym; j++) { 490c4762a1bSJed Brown for (i=info->xs; i<info->xs+info->xm; i++) { 491c4762a1bSJed Brown if (i == 0 || j == 0 || i == info->mx-1 || j == info->my-1) { 492c4762a1bSJed Brown lobj += PetscRealPart((hydhx + hxdhy)*x[j][i]*x[j][i]); 493c4762a1bSJed Brown } else { 494c4762a1bSJed Brown u = x[j][i]; 495c4762a1bSJed Brown uw = x[j][i-1]; 496c4762a1bSJed Brown ue = x[j][i+1]; 497c4762a1bSJed Brown un = x[j-1][i]; 498c4762a1bSJed Brown us = x[j+1][i]; 499c4762a1bSJed Brown 500c4762a1bSJed Brown if (i-1 == 0) uw = 0.; 501c4762a1bSJed Brown if (i+1 == info->mx-1) ue = 0.; 502c4762a1bSJed Brown if (j-1 == 0) un = 0.; 503c4762a1bSJed Brown if (j+1 == info->my-1) us = 0.; 504c4762a1bSJed Brown 505c4762a1bSJed Brown /* F[u] = 1/2\int_{\omega}\nabla^2u(x)*u(x)*dx */ 506c4762a1bSJed Brown 507c4762a1bSJed Brown uxux = u*(2.*u - ue - uw)*hydhx; 508c4762a1bSJed Brown uyuy = u*(2.*u - un - us)*hxdhy; 509c4762a1bSJed Brown 510c4762a1bSJed Brown lobj += PetscRealPart(0.5*(uxux + uyuy) - sc*PetscExpScalar(u)); 511c4762a1bSJed Brown } 512c4762a1bSJed Brown } 513c4762a1bSJed Brown } 514c4762a1bSJed Brown ierr = PetscLogFlops(12.0*info->ym*info->xm);CHKERRQ(ierr); 515c4762a1bSJed Brown ierr = MPI_Allreduce(&lobj,obj,1,MPIU_REAL,MPIU_SUM,comm);CHKERRQ(ierr); 516c4762a1bSJed Brown PetscFunctionReturn(0); 517c4762a1bSJed Brown } 518c4762a1bSJed Brown 519c4762a1bSJed Brown /* 520c4762a1bSJed Brown FormJacobianLocal - Evaluates Jacobian matrix on local process patch 521c4762a1bSJed Brown */ 522c4762a1bSJed Brown PetscErrorCode FormJacobianLocal(DMDALocalInfo *info,PetscScalar **x,Mat jac,Mat jacpre,AppCtx *user) 523c4762a1bSJed Brown { 524c4762a1bSJed Brown PetscErrorCode ierr; 525c4762a1bSJed Brown PetscInt i,j,k; 526c4762a1bSJed Brown MatStencil col[5],row; 527c4762a1bSJed Brown PetscScalar lambda,v[5],hx,hy,hxdhy,hydhx,sc; 528c4762a1bSJed Brown DM coordDA; 529c4762a1bSJed Brown Vec coordinates; 530c4762a1bSJed Brown DMDACoor2d **coords; 531c4762a1bSJed Brown 532c4762a1bSJed Brown PetscFunctionBeginUser; 533c4762a1bSJed Brown lambda = user->param; 534c4762a1bSJed Brown /* Extract coordinates */ 535c4762a1bSJed Brown ierr = DMGetCoordinateDM(info->da, &coordDA);CHKERRQ(ierr); 536c4762a1bSJed Brown ierr = DMGetCoordinates(info->da, &coordinates);CHKERRQ(ierr); 537c4762a1bSJed Brown ierr = DMDAVecGetArray(coordDA, coordinates, &coords);CHKERRQ(ierr); 538c4762a1bSJed Brown hx = info->xm > 1 ? PetscRealPart(coords[info->ys][info->xs+1].x) - PetscRealPart(coords[info->ys][info->xs].x) : 1.0; 539c4762a1bSJed Brown hy = info->ym > 1 ? PetscRealPart(coords[info->ys+1][info->xs].y) - PetscRealPart(coords[info->ys][info->xs].y) : 1.0; 540c4762a1bSJed Brown ierr = DMDAVecRestoreArray(coordDA, coordinates, &coords);CHKERRQ(ierr); 541c4762a1bSJed Brown hxdhy = hx/hy; 542c4762a1bSJed Brown hydhx = hy/hx; 543c4762a1bSJed Brown sc = hx*hy*lambda; 544c4762a1bSJed Brown 545c4762a1bSJed Brown 546c4762a1bSJed Brown /* 547c4762a1bSJed Brown Compute entries for the locally owned part of the Jacobian. 548c4762a1bSJed Brown - Currently, all PETSc parallel matrix formats are partitioned by 549c4762a1bSJed Brown contiguous chunks of rows across the processors. 550c4762a1bSJed Brown - Each processor needs to insert only elements that it owns 551c4762a1bSJed Brown locally (but any non-local elements will be sent to the 552c4762a1bSJed Brown appropriate processor during matrix assembly). 553c4762a1bSJed Brown - Here, we set all entries for a particular row at once. 554c4762a1bSJed Brown - We can set matrix entries either using either 555c4762a1bSJed Brown MatSetValuesLocal() or MatSetValues(), as discussed above. 556c4762a1bSJed Brown */ 557c4762a1bSJed Brown for (j=info->ys; j<info->ys+info->ym; j++) { 558c4762a1bSJed Brown for (i=info->xs; i<info->xs+info->xm; i++) { 559c4762a1bSJed Brown row.j = j; row.i = i; 560c4762a1bSJed Brown /* boundary points */ 561c4762a1bSJed Brown if (i == 0 || j == 0 || i == info->mx-1 || j == info->my-1) { 562c4762a1bSJed Brown v[0] = 2.0*(hydhx + hxdhy); 563c4762a1bSJed Brown ierr = MatSetValuesStencil(jacpre,1,&row,1,&row,v,INSERT_VALUES);CHKERRQ(ierr); 564c4762a1bSJed Brown } else { 565c4762a1bSJed Brown k = 0; 566c4762a1bSJed Brown /* interior grid points */ 567c4762a1bSJed Brown if (j-1 != 0) { 568c4762a1bSJed Brown v[k] = -hxdhy; 569c4762a1bSJed Brown col[k].j = j - 1; col[k].i = i; 570c4762a1bSJed Brown k++; 571c4762a1bSJed Brown } 572c4762a1bSJed Brown if (i-1 != 0) { 573c4762a1bSJed Brown v[k] = -hydhx; 574c4762a1bSJed Brown col[k].j = j; col[k].i = i-1; 575c4762a1bSJed Brown k++; 576c4762a1bSJed Brown } 577c4762a1bSJed Brown 578c4762a1bSJed Brown v[k] = 2.0*(hydhx + hxdhy) - sc*PetscExpScalar(x[j][i]); col[k].j = row.j; col[k].i = row.i; k++; 579c4762a1bSJed Brown 580c4762a1bSJed Brown if (i+1 != info->mx-1) { 581c4762a1bSJed Brown v[k] = -hydhx; 582c4762a1bSJed Brown col[k].j = j; col[k].i = i+1; 583c4762a1bSJed Brown k++; 584c4762a1bSJed Brown } 585c4762a1bSJed Brown if (j+1 != info->mx-1) { 586c4762a1bSJed Brown v[k] = -hxdhy; 587c4762a1bSJed Brown col[k].j = j + 1; col[k].i = i; 588c4762a1bSJed Brown k++; 589c4762a1bSJed Brown } 590c4762a1bSJed Brown ierr = MatSetValuesStencil(jacpre,1,&row,k,col,v,INSERT_VALUES);CHKERRQ(ierr); 591c4762a1bSJed Brown } 592c4762a1bSJed Brown } 593c4762a1bSJed Brown } 594c4762a1bSJed Brown 595c4762a1bSJed Brown /* 596c4762a1bSJed Brown Assemble matrix, using the 2-step process: 597c4762a1bSJed Brown MatAssemblyBegin(), MatAssemblyEnd(). 598c4762a1bSJed Brown */ 599c4762a1bSJed Brown ierr = MatAssemblyBegin(jacpre,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 600c4762a1bSJed Brown ierr = MatAssemblyEnd(jacpre,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 601c4762a1bSJed Brown /* 602c4762a1bSJed Brown Tell the matrix we will never add a new nonzero location to the 603c4762a1bSJed Brown matrix. If we do, it will generate an error. 604c4762a1bSJed Brown */ 605c4762a1bSJed Brown ierr = MatSetOption(jac,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr); 606c4762a1bSJed Brown PetscFunctionReturn(0); 607c4762a1bSJed Brown } 608c4762a1bSJed Brown 609c4762a1bSJed Brown PetscErrorCode FormFunctionMatlab(SNES snes,Vec X,Vec F,void *ptr) 610c4762a1bSJed Brown { 611*4e1ad211SJed Brown #if PetscDefined(HAVE_MATLAB_ENGINE) 612c4762a1bSJed Brown AppCtx *user = (AppCtx*)ptr; 613c4762a1bSJed Brown PetscErrorCode ierr; 614c4762a1bSJed Brown PetscInt Mx,My; 615c4762a1bSJed Brown PetscReal lambda,hx,hy; 616c4762a1bSJed Brown Vec localX,localF; 617c4762a1bSJed Brown MPI_Comm comm; 618c4762a1bSJed Brown DM da; 619c4762a1bSJed Brown 620c4762a1bSJed Brown PetscFunctionBeginUser; 621c4762a1bSJed Brown ierr = SNESGetDM(snes,&da);CHKERRQ(ierr); 622c4762a1bSJed Brown ierr = DMGetLocalVector(da,&localX);CHKERRQ(ierr); 623c4762a1bSJed Brown ierr = DMGetLocalVector(da,&localF);CHKERRQ(ierr); 624c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject)localX,"localX");CHKERRQ(ierr); 625c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject)localF,"localF");CHKERRQ(ierr); 626c4762a1bSJed Brown ierr = DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);CHKERRQ(ierr); 627c4762a1bSJed Brown 628c4762a1bSJed Brown lambda = user->param; 629c4762a1bSJed Brown hx = 1.0/(PetscReal)(Mx-1); 630c4762a1bSJed Brown hy = 1.0/(PetscReal)(My-1); 631c4762a1bSJed Brown 632c4762a1bSJed Brown ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr); 633c4762a1bSJed Brown /* 634c4762a1bSJed Brown Scatter ghost points to local vector,using the 2-step process 635c4762a1bSJed Brown DMGlobalToLocalBegin(),DMGlobalToLocalEnd(). 636c4762a1bSJed Brown By placing code between these two statements, computations can be 637c4762a1bSJed Brown done while messages are in transition. 638c4762a1bSJed Brown */ 639c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);CHKERRQ(ierr); 640c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);CHKERRQ(ierr); 641c4762a1bSJed Brown ierr = PetscMatlabEnginePut(PETSC_MATLAB_ENGINE_(comm),(PetscObject)localX);CHKERRQ(ierr); 642c4762a1bSJed Brown ierr = PetscMatlabEngineEvaluate(PETSC_MATLAB_ENGINE_(comm),"localF=ex5m(localX,%18.16e,%18.16e,%18.16e)",hx,hy,lambda);CHKERRQ(ierr); 643c4762a1bSJed Brown ierr = PetscMatlabEngineGet(PETSC_MATLAB_ENGINE_(comm),(PetscObject)localF);CHKERRQ(ierr); 644c4762a1bSJed Brown 645c4762a1bSJed Brown /* 646c4762a1bSJed Brown Insert values into global vector 647c4762a1bSJed Brown */ 648c4762a1bSJed Brown ierr = DMLocalToGlobalBegin(da,localF,INSERT_VALUES,F);CHKERRQ(ierr); 649c4762a1bSJed Brown ierr = DMLocalToGlobalEnd(da,localF,INSERT_VALUES,F);CHKERRQ(ierr); 650c4762a1bSJed Brown ierr = DMRestoreLocalVector(da,&localX);CHKERRQ(ierr); 651c4762a1bSJed Brown ierr = DMRestoreLocalVector(da,&localF);CHKERRQ(ierr); 652c4762a1bSJed Brown PetscFunctionReturn(0); 653*4e1ad211SJed Brown #else 654*4e1ad211SJed Brown return 0; /* Never called */ 655c4762a1bSJed Brown #endif 656*4e1ad211SJed Brown } 657c4762a1bSJed Brown 658c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 659c4762a1bSJed Brown /* 660c4762a1bSJed Brown Applies some sweeps on nonlinear Gauss-Seidel on each process 661c4762a1bSJed Brown 662c4762a1bSJed Brown */ 663c4762a1bSJed Brown PetscErrorCode NonlinearGS(SNES snes,Vec X, Vec B, void *ctx) 664c4762a1bSJed Brown { 665c4762a1bSJed Brown PetscInt i,j,k,Mx,My,xs,ys,xm,ym,its,tot_its,sweeps,l; 666c4762a1bSJed Brown PetscErrorCode ierr; 667c4762a1bSJed Brown PetscReal lambda,hx,hy,hxdhy,hydhx,sc; 668c4762a1bSJed Brown PetscScalar **x,**b,bij,F,F0=0,J,u,un,us,ue,eu,uw,uxx,uyy,y; 669c4762a1bSJed Brown PetscReal atol,rtol,stol; 670c4762a1bSJed Brown DM da; 671c4762a1bSJed Brown AppCtx *user; 672c4762a1bSJed Brown Vec localX,localB; 673c4762a1bSJed Brown 674c4762a1bSJed Brown PetscFunctionBeginUser; 675c4762a1bSJed Brown tot_its = 0; 676c4762a1bSJed Brown ierr = SNESNGSGetSweeps(snes,&sweeps);CHKERRQ(ierr); 677c4762a1bSJed Brown ierr = SNESNGSGetTolerances(snes,&atol,&rtol,&stol,&its);CHKERRQ(ierr); 678c4762a1bSJed Brown ierr = SNESGetDM(snes,&da);CHKERRQ(ierr); 679c4762a1bSJed Brown ierr = DMGetApplicationContext(da,(void**)&user);CHKERRQ(ierr); 680c4762a1bSJed Brown 681c4762a1bSJed Brown ierr = DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);CHKERRQ(ierr); 682c4762a1bSJed Brown 683c4762a1bSJed Brown lambda = user->param; 684c4762a1bSJed Brown hx = 1.0/(PetscReal)(Mx-1); 685c4762a1bSJed Brown hy = 1.0/(PetscReal)(My-1); 686c4762a1bSJed Brown sc = hx*hy*lambda; 687c4762a1bSJed Brown hxdhy = hx/hy; 688c4762a1bSJed Brown hydhx = hy/hx; 689c4762a1bSJed Brown 690c4762a1bSJed Brown 691c4762a1bSJed Brown ierr = DMGetLocalVector(da,&localX);CHKERRQ(ierr); 692c4762a1bSJed Brown if (B) { 693c4762a1bSJed Brown ierr = DMGetLocalVector(da,&localB);CHKERRQ(ierr); 694c4762a1bSJed Brown } 695c4762a1bSJed Brown for (l=0; l<sweeps; l++) { 696c4762a1bSJed Brown 697c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);CHKERRQ(ierr); 698c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);CHKERRQ(ierr); 699c4762a1bSJed Brown if (B) { 700c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(da,B,INSERT_VALUES,localB);CHKERRQ(ierr); 701c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(da,B,INSERT_VALUES,localB);CHKERRQ(ierr); 702c4762a1bSJed Brown } 703c4762a1bSJed Brown /* 704c4762a1bSJed Brown Get a pointer to vector data. 705c4762a1bSJed Brown - For default PETSc vectors, VecGetArray() returns a pointer to 706c4762a1bSJed Brown the data array. Otherwise, the routine is implementation dependent. 707c4762a1bSJed Brown - You MUST call VecRestoreArray() when you no longer need access to 708c4762a1bSJed Brown the array. 709c4762a1bSJed Brown */ 710c4762a1bSJed Brown ierr = DMDAVecGetArray(da,localX,&x);CHKERRQ(ierr); 711c4762a1bSJed Brown if (B) ierr = DMDAVecGetArray(da,localB,&b);CHKERRQ(ierr); 712c4762a1bSJed Brown /* 713c4762a1bSJed Brown Get local grid boundaries (for 2-dimensional DMDA): 714c4762a1bSJed Brown xs, ys - starting grid indices (no ghost points) 715c4762a1bSJed Brown xm, ym - widths of local grid (no ghost points) 716c4762a1bSJed Brown */ 717c4762a1bSJed Brown ierr = DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr); 718c4762a1bSJed Brown 719c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 720c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 721c4762a1bSJed Brown if (i == 0 || j == 0 || i == Mx-1 || j == My-1) { 722c4762a1bSJed Brown /* boundary conditions are all zero Dirichlet */ 723c4762a1bSJed Brown x[j][i] = 0.0; 724c4762a1bSJed Brown } else { 725c4762a1bSJed Brown if (B) bij = b[j][i]; 726c4762a1bSJed Brown else bij = 0.; 727c4762a1bSJed Brown 728c4762a1bSJed Brown u = x[j][i]; 729c4762a1bSJed Brown un = x[j-1][i]; 730c4762a1bSJed Brown us = x[j+1][i]; 731c4762a1bSJed Brown ue = x[j][i-1]; 732c4762a1bSJed Brown uw = x[j][i+1]; 733c4762a1bSJed Brown 734c4762a1bSJed Brown for (k=0; k<its; k++) { 735c4762a1bSJed Brown eu = PetscExpScalar(u); 736c4762a1bSJed Brown uxx = (2.0*u - ue - uw)*hydhx; 737c4762a1bSJed Brown uyy = (2.0*u - un - us)*hxdhy; 738c4762a1bSJed Brown F = uxx + uyy - sc*eu - bij; 739c4762a1bSJed Brown if (k == 0) F0 = F; 740c4762a1bSJed Brown J = 2.0*(hydhx + hxdhy) - sc*eu; 741c4762a1bSJed Brown y = F/J; 742c4762a1bSJed Brown u -= y; 743c4762a1bSJed Brown tot_its++; 744c4762a1bSJed Brown 745c4762a1bSJed Brown if (atol > PetscAbsReal(PetscRealPart(F)) || 746c4762a1bSJed Brown rtol*PetscAbsReal(PetscRealPart(F0)) > PetscAbsReal(PetscRealPart(F)) || 747c4762a1bSJed Brown stol*PetscAbsReal(PetscRealPart(u)) > PetscAbsReal(PetscRealPart(y))) { 748c4762a1bSJed Brown break; 749c4762a1bSJed Brown } 750c4762a1bSJed Brown } 751c4762a1bSJed Brown x[j][i] = u; 752c4762a1bSJed Brown } 753c4762a1bSJed Brown } 754c4762a1bSJed Brown } 755c4762a1bSJed Brown /* 756c4762a1bSJed Brown Restore vector 757c4762a1bSJed Brown */ 758c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da,localX,&x);CHKERRQ(ierr); 759c4762a1bSJed Brown ierr = DMLocalToGlobalBegin(da,localX,INSERT_VALUES,X);CHKERRQ(ierr); 760c4762a1bSJed Brown ierr = DMLocalToGlobalEnd(da,localX,INSERT_VALUES,X);CHKERRQ(ierr); 761c4762a1bSJed Brown } 762c4762a1bSJed Brown ierr = PetscLogFlops(tot_its*(21.0));CHKERRQ(ierr); 763c4762a1bSJed Brown ierr = DMRestoreLocalVector(da,&localX);CHKERRQ(ierr); 764c4762a1bSJed Brown if (B) { 765c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da,localB,&b);CHKERRQ(ierr); 766c4762a1bSJed Brown ierr = DMRestoreLocalVector(da,&localB);CHKERRQ(ierr); 767c4762a1bSJed Brown } 768c4762a1bSJed Brown PetscFunctionReturn(0); 769c4762a1bSJed Brown } 770c4762a1bSJed Brown 771c4762a1bSJed Brown /*TEST 772c4762a1bSJed Brown 773c4762a1bSJed Brown test: 774c4762a1bSJed Brown suffix: asm_0 775c4762a1bSJed Brown requires: !single 776c4762a1bSJed Brown args: -mms 1 -par 0.0 -snes_monitor_short -snes_converged_reason -snes_view -ksp_rtol 1.0e-9 -ksp_monitor_short -ksp_type richardson -pc_type asm -pc_asm_blocks 2 -pc_asm_overlap 0 -pc_asm_local_type additive -sub_pc_type lu 777c4762a1bSJed Brown 778c4762a1bSJed Brown test: 779c4762a1bSJed Brown suffix: msm_0 780c4762a1bSJed Brown requires: !single 781c4762a1bSJed Brown args: -mms 1 -par 0.0 -snes_monitor_short -snes_converged_reason -snes_view -ksp_rtol 1.0e-9 -ksp_monitor_short -ksp_type richardson -pc_type asm -pc_asm_blocks 2 -pc_asm_overlap 0 -pc_asm_local_type multiplicative -sub_pc_type lu 782c4762a1bSJed Brown 783c4762a1bSJed Brown test: 784c4762a1bSJed Brown suffix: asm_1 785c4762a1bSJed Brown requires: !single 786c4762a1bSJed Brown args: -mms 1 -par 0.0 -snes_monitor_short -snes_converged_reason -snes_view -ksp_rtol 1.0e-9 -ksp_monitor_short -ksp_type richardson -pc_type asm -pc_asm_blocks 2 -pc_asm_overlap 0 -pc_asm_local_type additive -sub_pc_type lu -da_grid_x 8 787c4762a1bSJed Brown 788c4762a1bSJed Brown test: 789c4762a1bSJed Brown suffix: msm_1 790c4762a1bSJed Brown requires: !single 791c4762a1bSJed Brown args: -mms 1 -par 0.0 -snes_monitor_short -snes_converged_reason -snes_view -ksp_rtol 1.0e-9 -ksp_monitor_short -ksp_type richardson -pc_type asm -pc_asm_blocks 2 -pc_asm_overlap 0 -pc_asm_local_type multiplicative -sub_pc_type lu -da_grid_x 8 792c4762a1bSJed Brown 793c4762a1bSJed Brown test: 794c4762a1bSJed Brown suffix: asm_2 795c4762a1bSJed Brown requires: !single 796c4762a1bSJed Brown args: -mms 1 -par 0.0 -snes_monitor_short -snes_converged_reason -snes_view -ksp_rtol 1.0e-9 -ksp_monitor_short -ksp_type richardson -pc_type asm -pc_asm_blocks 3 -pc_asm_overlap 0 -pc_asm_local_type additive -sub_pc_type lu -da_grid_x 8 797c4762a1bSJed Brown 798c4762a1bSJed Brown test: 799c4762a1bSJed Brown suffix: msm_2 800c4762a1bSJed Brown requires: !single 801c4762a1bSJed Brown args: -mms 1 -par 0.0 -snes_monitor_short -snes_converged_reason -snes_view -ksp_rtol 1.0e-9 -ksp_monitor_short -ksp_type richardson -pc_type asm -pc_asm_blocks 3 -pc_asm_overlap 0 -pc_asm_local_type multiplicative -sub_pc_type lu -da_grid_x 8 802c4762a1bSJed Brown 803c4762a1bSJed Brown test: 804c4762a1bSJed Brown suffix: asm_3 805c4762a1bSJed Brown requires: !single 806c4762a1bSJed Brown args: -mms 1 -par 0.0 -snes_monitor_short -snes_converged_reason -snes_view -ksp_rtol 1.0e-9 -ksp_monitor_short -ksp_type richardson -pc_type asm -pc_asm_blocks 4 -pc_asm_overlap 0 -pc_asm_local_type additive -sub_pc_type lu -da_grid_x 8 807c4762a1bSJed Brown 808c4762a1bSJed Brown test: 809c4762a1bSJed Brown suffix: msm_3 810c4762a1bSJed Brown requires: !single 811c4762a1bSJed Brown args: -mms 1 -par 0.0 -snes_monitor_short -snes_converged_reason -snes_view -ksp_rtol 1.0e-9 -ksp_monitor_short -ksp_type richardson -pc_type asm -pc_asm_blocks 4 -pc_asm_overlap 0 -pc_asm_local_type multiplicative -sub_pc_type lu -da_grid_x 8 812c4762a1bSJed Brown 813c4762a1bSJed Brown test: 814c4762a1bSJed Brown suffix: asm_4 815c4762a1bSJed Brown requires: !single 816c4762a1bSJed Brown nsize: 2 817c4762a1bSJed Brown args: -mms 1 -par 0.0 -snes_monitor_short -snes_converged_reason -snes_view -ksp_rtol 1.0e-9 -ksp_monitor_short -ksp_type richardson -pc_type asm -pc_asm_blocks 2 -pc_asm_overlap 0 -pc_asm_local_type additive -sub_pc_type lu -da_grid_x 8 818c4762a1bSJed Brown 819c4762a1bSJed Brown test: 820c4762a1bSJed Brown suffix: msm_4 821c4762a1bSJed Brown requires: !single 822c4762a1bSJed Brown nsize: 2 823c4762a1bSJed Brown args: -mms 1 -par 0.0 -snes_monitor_short -snes_converged_reason -snes_view -ksp_rtol 1.0e-9 -ksp_monitor_short -ksp_type richardson -pc_type asm -pc_asm_blocks 2 -pc_asm_overlap 0 -pc_asm_local_type multiplicative -sub_pc_type lu -da_grid_x 8 824c4762a1bSJed Brown 825c4762a1bSJed Brown test: 826c4762a1bSJed Brown suffix: asm_5 827c4762a1bSJed Brown requires: !single 828c4762a1bSJed Brown nsize: 2 829c4762a1bSJed Brown args: -mms 1 -par 0.0 -snes_monitor_short -snes_converged_reason -snes_view -ksp_rtol 1.0e-9 -ksp_monitor_short -ksp_type richardson -pc_type asm -pc_asm_blocks 4 -pc_asm_overlap 0 -pc_asm_local_type additive -sub_pc_type lu -da_grid_x 8 830c4762a1bSJed Brown 831c4762a1bSJed Brown test: 832c4762a1bSJed Brown suffix: msm_5 833c4762a1bSJed Brown requires: !single 834c4762a1bSJed Brown nsize: 2 835c4762a1bSJed Brown args: -mms 1 -par 0.0 -snes_monitor_short -snes_converged_reason -snes_view -ksp_rtol 1.0e-9 -ksp_monitor_short -ksp_type richardson -pc_type asm -pc_asm_blocks 4 -pc_asm_overlap 0 -pc_asm_local_type multiplicative -sub_pc_type lu -da_grid_x 8 836c4762a1bSJed Brown 837c4762a1bSJed Brown test: 838c4762a1bSJed Brown args: -snes_rtol 1.e-5 -pc_type mg -ksp_monitor_short -snes_view -pc_mg_levels 3 -pc_mg_galerkin pmat -da_grid_x 17 -da_grid_y 17 -mg_levels_ksp_monitor_short -mg_levels_ksp_norm_type unpreconditioned -snes_monitor_short -mg_levels_ksp_chebyshev_esteig 0.5,1.1 -mg_levels_pc_type sor -pc_mg_type full 839c4762a1bSJed Brown requires: !single 840c4762a1bSJed Brown 841c4762a1bSJed Brown test: 842c4762a1bSJed Brown suffix: 2 843c4762a1bSJed Brown args: -pc_type mg -ksp_converged_reason -snes_view -pc_mg_galerkin pmat -snes_grid_sequence 3 -mg_levels_ksp_norm_type unpreconditioned -snes_monitor_short -mg_levels_ksp_chebyshev_esteig 0.5,1.1 -mg_levels_pc_type sor -pc_mg_type full -ksp_atol -1. 844c4762a1bSJed Brown 845c4762a1bSJed Brown test: 846c4762a1bSJed Brown suffix: 3 847c4762a1bSJed Brown nsize: 2 848c4762a1bSJed Brown args: -snes_grid_sequence 2 -snes_mf_operator -snes_converged_reason -snes_view -pc_type mg -snes_atol -1 -snes_rtol 1.e-2 849c4762a1bSJed Brown filter: grep -v "otal number of function evaluations" 850c4762a1bSJed Brown 851c4762a1bSJed Brown test: 852c4762a1bSJed Brown suffix: 4 853c4762a1bSJed Brown nsize: 2 854c4762a1bSJed Brown args: -snes_grid_sequence 2 -snes_monitor_short -ksp_converged_reason -snes_converged_reason -snes_view -pc_type mg -snes_atol -1 -ksp_atol -1 855c4762a1bSJed Brown 856c4762a1bSJed Brown test: 857c4762a1bSJed Brown suffix: 5_anderson 858c4762a1bSJed Brown args: -da_grid_x 81 -da_grid_y 81 -snes_monitor_short -snes_max_it 50 -par 6.0 -snes_type anderson 859c4762a1bSJed Brown 860c4762a1bSJed Brown test: 861c4762a1bSJed Brown suffix: 5_aspin 862c4762a1bSJed Brown nsize: 4 863c4762a1bSJed Brown args: -snes_monitor_short -ksp_monitor_short -snes_converged_reason -da_refine 4 -da_overlap 3 -snes_type aspin -snes_view 864c4762a1bSJed Brown 865c4762a1bSJed Brown test: 866c4762a1bSJed Brown suffix: 5_broyden 867c4762a1bSJed Brown args: -da_grid_x 81 -da_grid_y 81 -snes_monitor_short -snes_max_it 50 -par 6.0 -snes_type qn -snes_qn_type broyden -snes_qn_m 10 868c4762a1bSJed Brown 869c4762a1bSJed Brown test: 870c4762a1bSJed Brown suffix: 5_fas 871c4762a1bSJed Brown args: -fas_coarse_snes_max_it 1 -fas_coarse_pc_type lu -fas_coarse_ksp_type preonly -snes_monitor_short -snes_type fas -fas_coarse_ksp_type richardson -da_refine 6 872c4762a1bSJed Brown requires: !single 873c4762a1bSJed Brown 874c4762a1bSJed Brown test: 875c4762a1bSJed Brown suffix: 5_fas_additive 876c4762a1bSJed Brown args: -fas_coarse_snes_max_it 1 -fas_coarse_pc_type lu -fas_coarse_ksp_type preonly -snes_monitor_short -snes_type fas -fas_coarse_ksp_type richardson -da_refine 6 -snes_fas_type additive -snes_max_it 50 877c4762a1bSJed Brown 878c4762a1bSJed Brown test: 879c4762a1bSJed Brown suffix: 5_fas_monitor 880c4762a1bSJed Brown args: -da_refine 1 -snes_type fas -snes_fas_monitor 881c4762a1bSJed Brown requires: !single 882c4762a1bSJed Brown 883c4762a1bSJed Brown test: 884c4762a1bSJed Brown suffix: 5_ls 885c4762a1bSJed Brown args: -da_grid_x 81 -da_grid_y 81 -snes_monitor_short -snes_max_it 50 -par 6.0 -snes_type newtonls 886c4762a1bSJed Brown 887c4762a1bSJed Brown test: 888c4762a1bSJed Brown suffix: 5_ls_sell_sor 889c4762a1bSJed Brown args: -da_grid_x 81 -da_grid_y 81 -snes_monitor_short -snes_max_it 50 -par 6.0 -snes_type newtonls -dm_mat_type sell -pc_type sor 890c4762a1bSJed Brown output_file: output/ex5_5_ls.out 891c4762a1bSJed Brown 892c4762a1bSJed Brown test: 893c4762a1bSJed Brown suffix: 5_nasm 894c4762a1bSJed Brown nsize: 4 895c4762a1bSJed Brown args: -snes_monitor_short -snes_converged_reason -da_refine 4 -da_overlap 3 -snes_type nasm -snes_nasm_type restrict -snes_max_it 10 896c4762a1bSJed Brown 897c4762a1bSJed Brown test: 898c4762a1bSJed Brown suffix: 5_ncg 899c4762a1bSJed Brown args: -da_grid_x 81 -da_grid_y 81 -snes_monitor_short -snes_max_it 50 -par 6.0 -snes_type ncg -snes_ncg_type fr 900c4762a1bSJed Brown 901c4762a1bSJed Brown test: 902c4762a1bSJed Brown suffix: 5_newton_asm_dmda 903c4762a1bSJed Brown nsize: 4 904c4762a1bSJed Brown args: -snes_monitor_short -ksp_monitor_short -snes_converged_reason -da_refine 4 -da_overlap 3 -snes_type newtonls -pc_type asm -pc_asm_dm_subdomains -malloc_dump 905c4762a1bSJed Brown requires: !single 906c4762a1bSJed Brown 907c4762a1bSJed Brown test: 908c4762a1bSJed Brown suffix: 5_newton_gasm_dmda 909c4762a1bSJed Brown nsize: 4 910c4762a1bSJed Brown args: -snes_monitor_short -ksp_monitor_short -snes_converged_reason -da_refine 4 -da_overlap 3 -snes_type newtonls -pc_type gasm -malloc_dump 911c4762a1bSJed Brown requires: !single 912c4762a1bSJed Brown 913c4762a1bSJed Brown test: 914c4762a1bSJed Brown suffix: 5_ngmres 915c4762a1bSJed Brown args: -da_grid_x 81 -da_grid_y 81 -snes_monitor_short -snes_max_it 50 -par 6.0 -snes_type ngmres -snes_ngmres_m 10 916c4762a1bSJed Brown 917c4762a1bSJed Brown test: 918c4762a1bSJed Brown suffix: 5_ngmres_fas 919c4762a1bSJed Brown args: -snes_rtol 1.e-4 -snes_type ngmres -npc_fas_coarse_snes_max_it 1 -npc_fas_coarse_snes_type newtonls -npc_fas_coarse_pc_type lu -npc_fas_coarse_ksp_type preonly -snes_ngmres_m 10 -snes_monitor_short -npc_snes_max_it 1 -npc_snes_type fas -npc_fas_coarse_ksp_type richardson -da_refine 6 920c4762a1bSJed Brown 921c4762a1bSJed Brown test: 922c4762a1bSJed Brown suffix: 5_ngmres_ngs 923c4762a1bSJed Brown args: -da_grid_x 81 -da_grid_y 81 -snes_monitor_short -snes_max_it 50 -par 6.0 -snes_type ngmres -npc_snes_type ngs -npc_snes_max_it 1 924c4762a1bSJed Brown 925c4762a1bSJed Brown test: 926c4762a1bSJed Brown suffix: 5_ngmres_nrichardson 927c4762a1bSJed Brown args: -da_grid_x 81 -da_grid_y 81 -snes_monitor_short -snes_max_it 50 -par 6.0 -snes_type ngmres -snes_ngmres_m 10 -npc_snes_type nrichardson -npc_snes_max_it 3 928c4762a1bSJed Brown output_file: output/ex5_5_ngmres_richardson.out 929c4762a1bSJed Brown 930c4762a1bSJed Brown test: 931c4762a1bSJed Brown suffix: 5_nrichardson 932c4762a1bSJed Brown args: -da_grid_x 81 -da_grid_y 81 -snes_monitor_short -snes_max_it 50 -par 6.0 -snes_type nrichardson 933c4762a1bSJed Brown 934c4762a1bSJed Brown test: 935c4762a1bSJed Brown suffix: 5_qn 936c4762a1bSJed Brown args: -da_grid_x 81 -da_grid_y 81 -snes_monitor_short -snes_max_it 50 -par 6.0 -snes_type qn -snes_linesearch_type cp -snes_qn_m 10 937c4762a1bSJed Brown 938c4762a1bSJed Brown test: 939c4762a1bSJed Brown suffix: 6 940c4762a1bSJed Brown nsize: 4 941c4762a1bSJed Brown args: -snes_converged_reason -ksp_converged_reason -da_grid_x 129 -da_grid_y 129 -pc_type mg -pc_mg_levels 8 -mg_levels_ksp_type chebyshev -mg_levels_ksp_chebyshev_esteig 0,0.5,0,1.1 -mg_levels_ksp_max_it 2 942c4762a1bSJed Brown 943c4762a1bSJed Brown test: 944c4762a1bSJed Brown requires: complex !single 945c4762a1bSJed Brown suffix: complex 946c4762a1bSJed Brown args: -snes_mf_operator -mat_mffd_complex -snes_monitor 947c4762a1bSJed Brown 948c4762a1bSJed Brown TEST*/ 949