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