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