1c4762a1bSJed Brown static char help[] = "Newton method to solve u'' + u^{2} = f, sequentially.\n\
2c4762a1bSJed Brown This example tests PCVPBJacobiSetBlocks().\n\n";
3c4762a1bSJed Brown
4c4762a1bSJed Brown /*
5c4762a1bSJed Brown Include "petscsnes.h" so that we can use SNES solvers. Note that this
6c4762a1bSJed Brown file automatically includes:
7c4762a1bSJed Brown petscsys.h - base PETSc routines petscvec.h - vectors
8c4762a1bSJed Brown petscmat.h - matrices
9c4762a1bSJed Brown petscis.h - index sets petscksp.h - Krylov subspace methods
10c4762a1bSJed Brown petscviewer.h - viewers petscpc.h - preconditioners
11c4762a1bSJed Brown petscksp.h - linear solvers
12c4762a1bSJed Brown */
13c4762a1bSJed Brown
14c4762a1bSJed Brown #include <petscsnes.h>
15c4762a1bSJed Brown
16c4762a1bSJed Brown /*
17c4762a1bSJed Brown User-defined routines
18c4762a1bSJed Brown */
19c4762a1bSJed Brown extern PetscErrorCode FormJacobian(SNES, Vec, Mat, Mat, void *);
20c4762a1bSJed Brown extern PetscErrorCode FormFunction(SNES, Vec, Vec, void *);
21c4762a1bSJed Brown extern PetscErrorCode FormInitialGuess(Vec);
22c4762a1bSJed Brown
main(int argc,char ** argv)23d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
24d71ae5a4SJacob Faibussowitsch {
25c4762a1bSJed Brown SNES snes; /* SNES context */
26c4762a1bSJed Brown Vec x, r, F, U; /* vectors */
27c4762a1bSJed Brown Mat J; /* Jacobian matrix */
28d7ff4d9fSStefano Zampini PetscInt its, n = 5, nb, maxit, maxf, *lens;
29c4762a1bSJed Brown PetscMPIInt size;
30c4762a1bSJed Brown PetscScalar h, xp, v, none = -1.0;
31c4762a1bSJed Brown PetscReal abstol, rtol, stol, norm;
32c4762a1bSJed Brown KSP ksp;
33c4762a1bSJed Brown PC pc;
34c4762a1bSJed Brown
35327415f7SBarry Smith PetscFunctionBeginUser;
36c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, NULL, help));
379566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
38d7ff4d9fSStefano Zampini PetscCheck(size == 1, PETSC_COMM_SELF, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only");
399566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &n, NULL));
40d7ff4d9fSStefano Zampini PetscCheck(n % 5 == 0, PETSC_COMM_SELF, PETSC_ERR_SUP, "n must be a multiple of 5");
41c4762a1bSJed Brown h = 1.0 / (n - 1);
42c4762a1bSJed Brown
43c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
44d7ff4d9fSStefano Zampini Create vector data structures
45c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
46c4762a1bSJed Brown /*
47c4762a1bSJed Brown Note that we form 1 vector from scratch and then duplicate as needed.
48c4762a1bSJed Brown */
499566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &x));
509566063dSJacob Faibussowitsch PetscCall(VecSetSizes(x, PETSC_DECIDE, n));
519566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(x));
529566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &r));
539566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &F));
549566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &U));
55c4762a1bSJed Brown
56c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
57d7ff4d9fSStefano Zampini Create matrix data structure
58c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
59c4762a1bSJed Brown
609566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &J));
619566063dSJacob Faibussowitsch PetscCall(MatSetSizes(J, PETSC_DECIDE, PETSC_DECIDE, n, n));
629566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(J));
639566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(J, 3, NULL));
64d7ff4d9fSStefano Zampini PetscCall(MatSetBlockSize(J, 5));
65d7ff4d9fSStefano Zampini
66d7ff4d9fSStefano Zampini nb = 3 * n / 5;
67d7ff4d9fSStefano Zampini PetscCall(PetscMalloc1(nb, &lens));
68d7ff4d9fSStefano Zampini for (PetscInt i = 0; i < nb / 3; i++) {
69d7ff4d9fSStefano Zampini lens[3 * i + 0] = 1;
70d7ff4d9fSStefano Zampini lens[3 * i + 1] = 2;
71d7ff4d9fSStefano Zampini lens[3 * i + 2] = 2;
72d7ff4d9fSStefano Zampini }
73d7ff4d9fSStefano Zampini PetscCall(MatSetVariableBlockSizes(J, nb, lens));
74d7ff4d9fSStefano Zampini PetscCall(PetscFree(lens));
75d7ff4d9fSStefano Zampini
76d7ff4d9fSStefano Zampini /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
77d7ff4d9fSStefano Zampini Create nonlinear solver context
78d7ff4d9fSStefano Zampini - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
79d7ff4d9fSStefano Zampini
80d7ff4d9fSStefano Zampini PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
81d7ff4d9fSStefano Zampini PetscCall(SNESGetKSP(snes, &ksp));
82d7ff4d9fSStefano Zampini PetscCall(KSPGetPC(ksp, &pc));
83d7ff4d9fSStefano Zampini PetscCall(PCSetType(pc, PCVPBJACOBI));
84d7ff4d9fSStefano Zampini
85d7ff4d9fSStefano Zampini /*
86d7ff4d9fSStefano Zampini Set function evaluation routine and vector
87d7ff4d9fSStefano Zampini */
88d7ff4d9fSStefano Zampini PetscCall(SNESSetFunction(snes, r, FormFunction, (void *)F));
89c4762a1bSJed Brown
90c4762a1bSJed Brown /*
91c4762a1bSJed Brown Set Jacobian matrix data structure and default Jacobian evaluation
92c4762a1bSJed Brown routine. User can override with:
93c4762a1bSJed Brown -snes_fd : default finite differencing approximation of Jacobian
94c4762a1bSJed Brown -snes_mf : matrix-free Newton-Krylov method with no preconditioning
95c4762a1bSJed Brown (unless user explicitly sets preconditioner)
967addb90fSBarry Smith -snes_mf_operator : form matrix used to construct the preconditioner as set by the user,
97c4762a1bSJed Brown but use matrix-free approx for Jacobian-vector
98c4762a1bSJed Brown products within Newton-Krylov method
99c4762a1bSJed Brown */
100c4762a1bSJed Brown
1019566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes, J, J, FormJacobian, NULL));
102c4762a1bSJed Brown
103c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
104c4762a1bSJed Brown Customize nonlinear solver; set runtime options
105c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
106c4762a1bSJed Brown
107c4762a1bSJed Brown /*
108c4762a1bSJed Brown Set names for some vectors to facilitate monitoring (optional)
109c4762a1bSJed Brown */
1109566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)x, "Approximate Solution"));
1119566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)U, "Exact Solution"));
112c4762a1bSJed Brown
113c4762a1bSJed Brown /*
114c4762a1bSJed Brown Set SNES/KSP/KSP/PC runtime options, e.g.,
115c4762a1bSJed Brown -snes_view -snes_monitor -ksp_type <ksp> -pc_type <pc>
116c4762a1bSJed Brown */
1179566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes));
118c4762a1bSJed Brown
119c4762a1bSJed Brown /*
120c4762a1bSJed Brown Print parameters used for convergence testing (optional) ... just
121c4762a1bSJed Brown to demonstrate this routine; this information is also printed with
122c4762a1bSJed Brown the option -snes_view
123c4762a1bSJed Brown */
1249566063dSJacob Faibussowitsch PetscCall(SNESGetTolerances(snes, &abstol, &rtol, &stol, &maxit, &maxf));
12563a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "atol=%g, rtol=%g, stol=%g, maxit=%" PetscInt_FMT ", maxf=%" PetscInt_FMT "\n", (double)abstol, (double)rtol, (double)stol, maxit, maxf));
126c4762a1bSJed Brown
127c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
128c4762a1bSJed Brown Initialize application:
129dd8e379bSPierre Jolivet Store right-hand side of PDE and exact solution
130c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
131c4762a1bSJed Brown
132c4762a1bSJed Brown xp = 0.0;
133d7ff4d9fSStefano Zampini for (PetscInt i = 0; i < n; i++) {
134c4762a1bSJed Brown v = 6.0 * xp + PetscPowScalar(xp + 1.e-12, 6.0); /* +1.e-12 is to prevent 0^6 */
1359566063dSJacob Faibussowitsch PetscCall(VecSetValues(F, 1, &i, &v, INSERT_VALUES));
136c4762a1bSJed Brown v = xp * xp * xp;
1379566063dSJacob Faibussowitsch PetscCall(VecSetValues(U, 1, &i, &v, INSERT_VALUES));
138c4762a1bSJed Brown xp += h;
139c4762a1bSJed Brown }
140c4762a1bSJed Brown
141c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
142c4762a1bSJed Brown Evaluate initial guess; then solve nonlinear system
143c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
144c4762a1bSJed Brown /*
145c4762a1bSJed Brown Note: The user should initialize the vector, x, with the initial guess
146c4762a1bSJed Brown for the nonlinear solver prior to calling SNESSolve(). In particular,
147c4762a1bSJed Brown to employ an initial guess of zero, the user should explicitly set
148c4762a1bSJed Brown this vector to zero by calling VecSet().
149c4762a1bSJed Brown */
1509566063dSJacob Faibussowitsch PetscCall(FormInitialGuess(x));
1519566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, NULL, x));
1529566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(snes, &its));
15363a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "number of SNES iterations = %" PetscInt_FMT "\n\n", its));
154c4762a1bSJed Brown
155c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
156c4762a1bSJed Brown Check solution and clean up
157c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
158c4762a1bSJed Brown
159c4762a1bSJed Brown /*
160c4762a1bSJed Brown Check the error
161c4762a1bSJed Brown */
1629566063dSJacob Faibussowitsch PetscCall(VecAXPY(x, none, U));
1639566063dSJacob Faibussowitsch PetscCall(VecNorm(x, NORM_2, &norm));
16463a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Norm of error %g, Iterations %" PetscInt_FMT "\n", (double)norm, its));
165c4762a1bSJed Brown
166c4762a1bSJed Brown /*
167c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they
168c4762a1bSJed Brown are no longer needed.
169c4762a1bSJed Brown */
1709371c9d4SSatish Balay PetscCall(VecDestroy(&x));
1719371c9d4SSatish Balay PetscCall(VecDestroy(&r));
1729371c9d4SSatish Balay PetscCall(VecDestroy(&U));
1739371c9d4SSatish Balay PetscCall(VecDestroy(&F));
1749371c9d4SSatish Balay PetscCall(MatDestroy(&J));
1759371c9d4SSatish Balay PetscCall(SNESDestroy(&snes));
1769566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
177b122ec5aSJacob Faibussowitsch return 0;
178c4762a1bSJed Brown }
179f6dfbefdSBarry Smith
180c4762a1bSJed Brown /*
181c4762a1bSJed Brown FormInitialGuess - Computes initial guess.
182c4762a1bSJed Brown
183c4762a1bSJed Brown Input/Output Parameter:
184c4762a1bSJed Brown . x - the solution vector
185c4762a1bSJed Brown */
FormInitialGuess(Vec x)186d71ae5a4SJacob Faibussowitsch PetscErrorCode FormInitialGuess(Vec x)
187d71ae5a4SJacob Faibussowitsch {
188c4762a1bSJed Brown PetscScalar pfive = .50;
1893ba16761SJacob Faibussowitsch
1903ba16761SJacob Faibussowitsch PetscFunctionBeginUser;
1919566063dSJacob Faibussowitsch PetscCall(VecSet(x, pfive));
1923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
193c4762a1bSJed Brown }
194f6dfbefdSBarry Smith
195c4762a1bSJed Brown /*
196c4762a1bSJed Brown FormFunction - Evaluates nonlinear function, F(x).
197c4762a1bSJed Brown
198c4762a1bSJed Brown Input Parameters:
199c4762a1bSJed Brown . snes - the SNES context
200c4762a1bSJed Brown . x - input vector
201c4762a1bSJed Brown . ctx - optional user-defined context, as set by SNESSetFunction()
202c4762a1bSJed Brown
203c4762a1bSJed Brown Output Parameter:
204c4762a1bSJed Brown . f - function vector
205c4762a1bSJed Brown
206c4762a1bSJed Brown Note:
207c4762a1bSJed Brown The user-defined context can contain any application-specific data
208c4762a1bSJed Brown needed for the function evaluation (such as various parameters, work
209c4762a1bSJed Brown vectors, and grid information). In this program the context is just
210dd8e379bSPierre Jolivet a vector containing the right-hand side of the discretized PDE.
211c4762a1bSJed Brown */
212c4762a1bSJed Brown
FormFunction(SNES snes,Vec x,Vec f,PetscCtx ctx)213*2a8381b2SBarry Smith PetscErrorCode FormFunction(SNES snes, Vec x, Vec f, PetscCtx ctx)
214d71ae5a4SJacob Faibussowitsch {
215c4762a1bSJed Brown Vec g = (Vec)ctx;
216c4762a1bSJed Brown const PetscScalar *xx, *gg;
217c4762a1bSJed Brown PetscScalar *ff, d;
218c4762a1bSJed Brown PetscInt i, n;
219c4762a1bSJed Brown
2203ba16761SJacob Faibussowitsch PetscFunctionBeginUser;
221c4762a1bSJed Brown /*
222c4762a1bSJed Brown Get pointers to vector data.
223c4762a1bSJed Brown - For default PETSc vectors, VecGetArray() returns a pointer to
224c4762a1bSJed Brown the data array. Otherwise, the routine is implementation dependent.
225c4762a1bSJed Brown - You MUST call VecRestoreArray() when you no longer need access to
226c4762a1bSJed Brown the array.
227c4762a1bSJed Brown */
2289566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &xx));
2299566063dSJacob Faibussowitsch PetscCall(VecGetArray(f, &ff));
2309566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(g, &gg));
231c4762a1bSJed Brown
232c4762a1bSJed Brown /*
233c4762a1bSJed Brown Compute function
234c4762a1bSJed Brown */
2359566063dSJacob Faibussowitsch PetscCall(VecGetSize(x, &n));
2369371c9d4SSatish Balay d = (PetscReal)(n - 1);
2379371c9d4SSatish Balay d = d * d;
238c4762a1bSJed Brown ff[0] = xx[0];
239c4762a1bSJed Brown for (i = 1; i < n - 1; i++) ff[i] = d * (xx[i - 1] - 2.0 * xx[i] + xx[i + 1]) + xx[i] * xx[i] - gg[i];
240c4762a1bSJed Brown ff[n - 1] = xx[n - 1] - 1.0;
241c4762a1bSJed Brown
242c4762a1bSJed Brown /*
243c4762a1bSJed Brown Restore vectors
244c4762a1bSJed Brown */
2459566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &xx));
2469566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(f, &ff));
2479566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(g, &gg));
2483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
249c4762a1bSJed Brown }
250f6dfbefdSBarry Smith
251c4762a1bSJed Brown /*
252c4762a1bSJed Brown FormJacobian - Evaluates Jacobian matrix.
253c4762a1bSJed Brown
254c4762a1bSJed Brown Input Parameters:
255c4762a1bSJed Brown . snes - the SNES context
256c4762a1bSJed Brown . x - input vector
257c4762a1bSJed Brown . dummy - optional user-defined context (not used here)
258c4762a1bSJed Brown
259c4762a1bSJed Brown Output Parameters:
260c4762a1bSJed Brown . jac - Jacobian matrix
2617addb90fSBarry Smith . B - optionally different matrix used to construct the preconditioner
262c4762a1bSJed Brown
263c4762a1bSJed Brown */
264c4762a1bSJed Brown
FormJacobian(SNES snes,Vec x,Mat jac,Mat B,void * dummy)265d71ae5a4SJacob Faibussowitsch PetscErrorCode FormJacobian(SNES snes, Vec x, Mat jac, Mat B, void *dummy)
266d71ae5a4SJacob Faibussowitsch {
267c4762a1bSJed Brown const PetscScalar *xx;
268c4762a1bSJed Brown PetscScalar A[3], d;
269c4762a1bSJed Brown PetscInt i, n, j[3];
270c4762a1bSJed Brown
2713ba16761SJacob Faibussowitsch PetscFunctionBeginUser;
272c4762a1bSJed Brown /*
273c4762a1bSJed Brown Get pointer to vector data
274c4762a1bSJed Brown */
2759566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &xx));
276c4762a1bSJed Brown
277c4762a1bSJed Brown /*
278c4762a1bSJed Brown Compute Jacobian entries and insert into matrix.
279c4762a1bSJed Brown - Note that in this case we set all elements for a particular
280c4762a1bSJed Brown row at once.
281c4762a1bSJed Brown */
2829566063dSJacob Faibussowitsch PetscCall(VecGetSize(x, &n));
2839371c9d4SSatish Balay d = (PetscReal)(n - 1);
2849371c9d4SSatish Balay d = d * d;
285c4762a1bSJed Brown
286c4762a1bSJed Brown /*
287c4762a1bSJed Brown Interior grid points
288c4762a1bSJed Brown */
289c4762a1bSJed Brown for (i = 1; i < n - 1; i++) {
2909371c9d4SSatish Balay j[0] = i - 1;
2919371c9d4SSatish Balay j[1] = i;
2929371c9d4SSatish Balay j[2] = i + 1;
2938d47944aSStefano Zampini A[0] = d;
2949371c9d4SSatish Balay A[1] = -2.0 * d + 2.0 * xx[i];
2958d47944aSStefano Zampini A[2] = d;
2969566063dSJacob Faibussowitsch PetscCall(MatSetValues(B, 1, &i, 3, j, A, INSERT_VALUES));
297c4762a1bSJed Brown }
298c4762a1bSJed Brown
299c4762a1bSJed Brown /*
300c4762a1bSJed Brown Boundary points
301c4762a1bSJed Brown */
3029371c9d4SSatish Balay i = 0;
3039371c9d4SSatish Balay A[0] = 1.0;
304c4762a1bSJed Brown
3059566063dSJacob Faibussowitsch PetscCall(MatSetValues(B, 1, &i, 1, &i, A, INSERT_VALUES));
306c4762a1bSJed Brown
3079371c9d4SSatish Balay i = n - 1;
3089371c9d4SSatish Balay A[0] = 1.0;
309c4762a1bSJed Brown
3109566063dSJacob Faibussowitsch PetscCall(MatSetValues(B, 1, &i, 1, &i, A, INSERT_VALUES));
311c4762a1bSJed Brown
312c4762a1bSJed Brown /*
313c4762a1bSJed Brown Restore vector
314c4762a1bSJed Brown */
3159566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &xx));
316c4762a1bSJed Brown
317c4762a1bSJed Brown /*
318c4762a1bSJed Brown Assemble matrix
319c4762a1bSJed Brown */
3209566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
3219566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
322c4762a1bSJed Brown if (jac != B) {
3239566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY));
3249566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY));
325c4762a1bSJed Brown }
3263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
327c4762a1bSJed Brown }
328c4762a1bSJed Brown
329c4762a1bSJed Brown /*TEST
330c4762a1bSJed Brown
331f1be3500SJunchao Zhang testset:
332c4762a1bSJed Brown args: -snes_monitor_short -snes_view -ksp_monitor
333f1be3500SJunchao Zhang output_file: output/ex5_1.out
334f1be3500SJunchao Zhang filter: grep -v "type: seqaij"
335f1be3500SJunchao Zhang
336f1be3500SJunchao Zhang test:
337f1be3500SJunchao Zhang suffix: 1
338f1be3500SJunchao Zhang
339f1be3500SJunchao Zhang test:
340f1be3500SJunchao Zhang suffix: cuda
341f1be3500SJunchao Zhang requires: cuda
342f1be3500SJunchao Zhang args: -mat_type aijcusparse -vec_type cuda
343f1be3500SJunchao Zhang
344f1be3500SJunchao Zhang test:
345f1be3500SJunchao Zhang suffix: kok
346f1be3500SJunchao Zhang requires: kokkos_kernels
347f1be3500SJunchao Zhang args: -mat_type aijkokkos -vec_type kokkos
348c4762a1bSJed Brown
34993d70b8aSPierre Jolivet # this is just a test for SNESKSPTRANSPOSEONLY and KSPSolveTranspose() to behave properly
350c4762a1bSJed Brown # the solution is wrong on purpose
351c4762a1bSJed Brown test:
352c4762a1bSJed Brown requires: !single !complex
353c4762a1bSJed Brown suffix: transpose_only
354c4762a1bSJed Brown args: -snes_monitor_short -snes_view -ksp_monitor -snes_type ksptransposeonly -pc_type ilu -snes_test_jacobian -snes_test_jacobian_view -ksp_view_rhs -ksp_view_solution -ksp_view_mat_explicit -ksp_view_preconditioned_operator_explicit
355c4762a1bSJed Brown
35693d70b8aSPierre Jolivet test:
35793d70b8aSPierre Jolivet requires: mumps
35893d70b8aSPierre Jolivet suffix: mumps
35993d70b8aSPierre Jolivet args: -pc_type lu -pc_factor_mat_solver_type mumps -mat_mumps_icntl_15 1 -snes_monitor_short -ksp_monitor
36093d70b8aSPierre Jolivet
361d7ff4d9fSStefano Zampini test:
362d7ff4d9fSStefano Zampini suffix: fieldsplit_1
363d7ff4d9fSStefano Zampini args: -snes_monitor_short -snes_view -ksp_monitor -pc_type fieldsplit -pc_fieldsplit_0_fields 0,1 -pc_fieldsplit_1_fields 2,3,4
364d7ff4d9fSStefano Zampini
365d7ff4d9fSStefano Zampini test:
366d7ff4d9fSStefano Zampini suffix: fieldsplit_2
367d7ff4d9fSStefano Zampini args: -snes_monitor_short -snes_view -ksp_monitor -pc_type fieldsplit -pc_fieldsplit_0_fields 1,0 -pc_fieldsplit_1_fields 2,3,4
368d7ff4d9fSStefano Zampini
369c4762a1bSJed Brown TEST*/
370