Actual source code: ex40.c
slepc-main 2024-11-09
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
7: SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9: */
11: static char help[] = "Checking the definite property in quadratic symmetric eigenproblem.\n\n"
12: "The command line options are:\n"
13: " -n <n> ... dimension of the matrices.\n"
14: " -transform... whether to transform to a hyperbolic problem or not.\n"
15: " -nonhyperbolic... to test with a modified (definite) problem that is not hyperbolic.\n\n";
17: #include <slepcpep.h>
19: /*
20: This example is based on spring.c, for fixed values mu=1,tau=10,kappa=5
22: The transformations are based on the method proposed in [Niendorf and Voss, LAA 2010].
23: */
25: PetscErrorCode QEPDefiniteTransformGetMatrices(PEP,PetscBool,PetscReal,PetscReal,Mat[3]);
26: PetscErrorCode QEPDefiniteTransformMap(PetscBool,PetscReal,PetscReal,PetscInt,PetscScalar*,PetscBool);
27: PetscErrorCode QEPDefiniteCheckError(Mat*,PEP,PetscBool,PetscReal,PetscReal);
28: PetscErrorCode TransformMatricesMoebius(Mat[3],MatStructure,PetscReal,PetscReal,PetscReal,PetscReal,Mat[3]);
30: int main(int argc,char **argv)
31: {
32: Mat M,C,K,*Op,A[3],At[3],B[3]; /* problem matrices */
33: PEP pep; /* polynomial eigenproblem solver context */
34: ST st; /* spectral transformation context */
35: KSP ksp;
36: PC pc;
37: PEPProblemType type;
38: PetscBool terse,transform=PETSC_FALSE,nohyp=PETSC_FALSE;
39: PetscInt n=100,Istart,Iend,i,def=0,hyp;
40: PetscReal muu=1,tau=10,kappa=5,inta,intb;
41: PetscReal alpha,beta,xi,mu,at[2]={0.0,0.0},c=.857,s;
42: PetscScalar target,targett,ats[2];
44: PetscFunctionBeginUser;
45: PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
47: PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
48: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nPEP example that checks definite property, n=%" PetscInt_FMT "\n\n",n));
50: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
51: Compute the matrices that define the eigensystem, (k^2*M+k*C+K)x=0
52: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
54: /* K is a tridiagonal */
55: PetscCall(MatCreate(PETSC_COMM_WORLD,&K));
56: PetscCall(MatSetSizes(K,PETSC_DECIDE,PETSC_DECIDE,n,n));
57: PetscCall(MatSetFromOptions(K));
59: PetscCall(MatGetOwnershipRange(K,&Istart,&Iend));
60: for (i=Istart;i<Iend;i++) {
61: if (i>0) PetscCall(MatSetValue(K,i,i-1,-kappa,INSERT_VALUES));
62: PetscCall(MatSetValue(K,i,i,kappa*3.0,INSERT_VALUES));
63: if (i<n-1) PetscCall(MatSetValue(K,i,i+1,-kappa,INSERT_VALUES));
64: }
66: PetscCall(MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY));
67: PetscCall(MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY));
69: /* C is a tridiagonal */
70: PetscCall(MatCreate(PETSC_COMM_WORLD,&C));
71: PetscCall(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,n,n));
72: PetscCall(MatSetFromOptions(C));
74: PetscCall(MatGetOwnershipRange(C,&Istart,&Iend));
75: for (i=Istart;i<Iend;i++) {
76: if (i>0) PetscCall(MatSetValue(C,i,i-1,-tau,INSERT_VALUES));
77: PetscCall(MatSetValue(C,i,i,tau*3.0,INSERT_VALUES));
78: if (i<n-1) PetscCall(MatSetValue(C,i,i+1,-tau,INSERT_VALUES));
79: }
81: PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
82: PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
84: /* M is a diagonal matrix */
85: PetscCall(MatCreate(PETSC_COMM_WORLD,&M));
86: PetscCall(MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,n,n));
87: PetscCall(MatSetFromOptions(M));
88: PetscCall(MatGetOwnershipRange(M,&Istart,&Iend));
89: for (i=Istart;i<Iend;i++) PetscCall(MatSetValue(M,i,i,muu,INSERT_VALUES));
90: PetscCall(MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY));
91: PetscCall(MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY));
93: PetscCall(PetscOptionsGetBool(NULL,NULL,"-nonhyperbolic",&nohyp,NULL));
94: A[0] = K; A[1] = C; A[2] = M;
95: if (nohyp) {
96: s = c*.6;
97: PetscCall(TransformMatricesMoebius(A,UNKNOWN_NONZERO_PATTERN,c,s,-s,c,At));
98: for (i=0;i<3;i++) PetscCall(MatDestroy(&A[i]));
99: Op = At;
100: } else Op = A;
102: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
103: Create the eigensolver and solve the problem
104: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
106: /*
107: Create eigensolver context
108: */
109: PetscCall(PEPCreate(PETSC_COMM_WORLD,&pep));
110: PetscCall(PEPSetProblemType(pep,PEP_HERMITIAN));
111: PetscCall(PEPSetType(pep,PEPSTOAR));
112: /*
113: Set operators and set problem type
114: */
115: PetscCall(PEPSetOperators(pep,3,Op));
117: /*
118: Set shift-and-invert with Cholesky; select MUMPS if available
119: */
120: PetscCall(PEPGetST(pep,&st));
121: PetscCall(STGetKSP(st,&ksp));
122: PetscCall(KSPSetType(ksp,KSPPREONLY));
123: PetscCall(KSPGetPC(ksp,&pc));
124: PetscCall(PCSetType(pc,PCCHOLESKY));
126: /*
127: Use MUMPS if available.
128: Note that in complex scalars we cannot use MUMPS for spectrum slicing,
129: because MatGetInertia() is not available in that case.
130: */
131: #if defined(PETSC_HAVE_MUMPS) && !defined(PETSC_USE_COMPLEX)
132: PetscCall(PCFactorSetMatSolverType(pc,MATSOLVERMUMPS));
133: /*
134: Add several MUMPS options (see ex43.c for a better way of setting them in program):
135: '-st_mat_mumps_icntl_13 1': turn off ScaLAPACK for matrix inertia
136: */
137: PetscCall(PetscOptionsInsertString(NULL,"-st_mat_mumps_icntl_13 1 -st_mat_mumps_icntl_24 1 -st_mat_mumps_cntl_3 1e-12"));
138: #endif
140: /*
141: Set solver parameters at runtime
142: */
143: PetscCall(PEPSetFromOptions(pep));
145: PetscCall(PetscOptionsGetBool(NULL,NULL,"-transform",&transform,NULL));
146: if (transform) {
147: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
148: Check if the problem is definite
149: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
150: PetscCall(PEPCheckDefiniteQEP(pep,&xi,&mu,&def,&hyp));
151: switch (def) {
152: case 1:
153: if (hyp==1) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Hyperbolic Problem xi=%g\n",(double)xi));
154: else PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Definite Problem xi=%g mu=%g\n",(double)xi,(double)mu));
155: break;
156: case -1:
157: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Not Definite Problem\n"));
158: break;
159: default:
160: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Cannot determine definiteness\n"));
161: break;
162: }
164: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
165: Transform the QEP to have a definite inner product in the linearization
166: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
167: if (def==1) {
168: PetscCall(QEPDefiniteTransformGetMatrices(pep,hyp==1?PETSC_TRUE:PETSC_FALSE,xi,mu,B));
169: PetscCall(PEPSetOperators(pep,3,B));
170: PetscCall(PEPGetTarget(pep,&target));
171: targett = target;
172: PetscCall(QEPDefiniteTransformMap(hyp==1?PETSC_TRUE:PETSC_FALSE,xi,mu,1,&targett,PETSC_FALSE));
173: PetscCall(PEPSetTarget(pep,targett));
174: PetscCall(PEPGetProblemType(pep,&type));
175: PetscCall(PEPSetProblemType(pep,PEP_HYPERBOLIC));
176: PetscCall(PEPSTOARGetLinearization(pep,&alpha,&beta));
177: PetscCall(PEPSTOARSetLinearization(pep,1.0,0.0));
178: PetscCall(PEPGetInterval(pep,&inta,&intb));
179: if (inta!=intb) {
180: ats[0] = inta; ats[1] = intb;
181: PetscCall(QEPDefiniteTransformMap(hyp==1?PETSC_TRUE:PETSC_FALSE,xi,mu,2,ats,PETSC_FALSE));
182: at[0] = PetscRealPart(ats[0]); at[1] = PetscRealPart(ats[1]);
183: if (at[0]<at[1]) PetscCall(PEPSetInterval(pep,at[0],at[1]));
184: else PetscCall(PEPSetInterval(pep,PETSC_MIN_REAL,at[1]));
185: }
186: }
187: }
189: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
190: Solve the eigensystem
191: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
192: PetscCall(PEPSolve(pep));
194: /* show detailed info unless -terse option is given by user */
195: if (def!=1) {
196: PetscCall(PetscOptionsHasName(NULL,NULL,"-terse",&terse));
197: if (terse) PetscCall(PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL));
198: else {
199: PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL));
200: PetscCall(PEPConvergedReasonView(pep,PETSC_VIEWER_STDOUT_WORLD));
201: PetscCall(PEPErrorView(pep,PEP_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD));
202: PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
203: }
204: } else {
205: /* Map the solution */
206: PetscCall(PEPConvergedReasonView(pep,PETSC_VIEWER_STDOUT_WORLD));
207: PetscCall(QEPDefiniteCheckError(Op,pep,hyp==1?PETSC_TRUE:PETSC_FALSE,xi,mu));
208: for (i=0;i<3;i++) PetscCall(MatDestroy(B+i));
209: }
210: if (at[0]>at[1]) {
211: PetscCall(PEPSetInterval(pep,at[0],PETSC_MAX_REAL));
212: PetscCall(PEPSolve(pep));
213: PetscCall(PEPConvergedReasonView(pep,PETSC_VIEWER_STDOUT_WORLD));
214: /* Map the solution */
215: PetscCall(QEPDefiniteCheckError(Op,pep,hyp==1?PETSC_TRUE:PETSC_FALSE,xi,mu));
216: }
217: if (def==1) {
218: PetscCall(PEPSetTarget(pep,target));
219: PetscCall(PEPSetProblemType(pep,type));
220: PetscCall(PEPSTOARSetLinearization(pep,alpha,beta));
221: if (inta!=intb) PetscCall(PEPSetInterval(pep,inta,intb));
222: }
224: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
225: Clean up
226: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
227: PetscCall(PEPDestroy(&pep));
228: for (i=0;i<3;i++) PetscCall(MatDestroy(Op+i));
229: PetscCall(SlepcFinalize());
230: return 0;
231: }
233: /* ------------------------------------------------------------------- */
234: /*
235: QEPDefiniteTransformMap_Initial - map a scalar value with a certain Moebius transform
237: a theta + b
238: lambda = --------------
239: c theta + d
241: Input:
242: xi,mu: real values such that Q(xi)<0 and Q(mu)>0
243: hyperbolic: if true the problem is assumed hyperbolic (mu is not used)
244: Input/Output:
245: val (array of length n)
246: if backtransform=true returns lambda from theta, else returns theta from lambda
247: */
248: static PetscErrorCode QEPDefiniteTransformMap_Initial(PetscBool hyperbolic,PetscReal xi,PetscReal mu,PetscInt n,PetscScalar *val,PetscBool backtransform)
249: {
250: PetscInt i;
251: PetscReal a,b,c,d,s;
253: PetscFunctionBegin;
254: if (hyperbolic) { a = 1.0; b = xi; c =0.0; d = 1.0; }
255: else { a = mu; b = mu*xi-1; c = 1.0; d = xi+mu; }
256: if (!backtransform) { s = a; a = -d; d = -s; }
257: for (i=0;i<n;i++) {
258: if (PetscRealPart(val[i]) >= PETSC_MAX_REAL || PetscRealPart(val[i]) <= PETSC_MIN_REAL) val[i] = a/c;
259: else if (val[i] == -d/c) val[i] = PETSC_MAX_REAL;
260: else val[i] = (a*val[i]+b)/(c*val[i]+d);
261: }
262: PetscFunctionReturn(PETSC_SUCCESS);
263: }
265: /* ------------------------------------------------------------------- */
266: /*
267: QEPDefiniteTransformMap - perform the mapping if the problem is hyperbolic, otherwise
268: modify the value of xi in advance
269: */
270: PetscErrorCode QEPDefiniteTransformMap(PetscBool hyperbolic,PetscReal xi,PetscReal mu,PetscInt n,PetscScalar *val,PetscBool backtransform)
271: {
272: PetscReal xit;
273: PetscScalar alpha;
275: PetscFunctionBegin;
276: xit = xi;
277: if (!hyperbolic) {
278: alpha = xi;
279: PetscCall(QEPDefiniteTransformMap_Initial(PETSC_FALSE,0.0,mu,1,&alpha,PETSC_FALSE));
280: xit = PetscRealPart(alpha);
281: }
282: PetscCall(QEPDefiniteTransformMap_Initial(hyperbolic,xit,mu,n,val,backtransform));
283: PetscFunctionReturn(PETSC_SUCCESS);
284: }
286: /* ------------------------------------------------------------------- */
287: /*
288: TransformMatricesMoebius - transform the coefficient matrices of a QEP
290: Input:
291: A: coefficient matrices of the original QEP
292: a,b,c,d: parameters of the Moebius transform
293: str: structure flag for MatAXPY operations
294: Output:
295: B: transformed matrices
296: */
297: PetscErrorCode TransformMatricesMoebius(Mat A[3],MatStructure str,PetscReal a,PetscReal b,PetscReal c,PetscReal d,Mat B[3])
298: {
299: PetscInt i,k;
300: PetscReal cf[9];
302: PetscFunctionBegin;
303: for (i=0;i<3;i++) PetscCall(MatDuplicate(A[2],MAT_COPY_VALUES,&B[i]));
304: /* Ct = b*b*A+b*d*B+d*d*C */
305: cf[0] = d*d; cf[1] = b*d; cf[2] = b*b;
306: /* Bt = 2*a*b*A+(b*c+a*d)*B+2*c*d*C*/
307: cf[3] = 2*c*d; cf[4] = b*c+a*d; cf[5] = 2*a*b;
308: /* At = a*a*A+a*c*B+c*c*C */
309: cf[6] = c*c; cf[7] = a*c; cf[8] = a*a;
310: for (k=0;k<3;k++) {
311: PetscCall(MatScale(B[k],cf[k*3+2]));
312: for (i=0;i<2;i++) PetscCall(MatAXPY(B[k],cf[3*k+i],A[i],str));
313: }
314: PetscFunctionReturn(PETSC_SUCCESS);
315: }
317: /* ------------------------------------------------------------------- */
318: /*
319: QEPDefiniteTransformGetMatrices - given a PEP of degree 2, transform the three
320: matrices with TransformMatricesMoebius
322: Input:
323: pep: polynomial eigenproblem to be transformed, with Q(.) being the quadratic polynomial
324: xi,mu: real values such that Q(xi)<0 and Q(mu)>0
325: hyperbolic: if true the problem is assumed hyperbolic (mu is not used)
326: Output:
327: T: coefficient matrices of the transformed polynomial
328: */
329: PetscErrorCode QEPDefiniteTransformGetMatrices(PEP pep,PetscBool hyperbolic,PetscReal xi,PetscReal mu,Mat T[3])
330: {
331: MatStructure str;
332: ST st;
333: PetscInt i;
334: PetscReal a,b,c,d;
335: PetscScalar xit;
336: Mat A[3];
338: PetscFunctionBegin;
339: for (i=2;i>=0;i--) PetscCall(PEPGetOperators(pep,i,&A[i]));
340: if (hyperbolic) { a = 1.0; b = xi; c =0.0; d = 1.0; }
341: else {
342: xit = xi;
343: PetscCall(QEPDefiniteTransformMap_Initial(PETSC_FALSE,0.0,mu,1,&xit,PETSC_FALSE));
344: a = mu; b = mu*PetscRealPart(xit)-1.0; c = 1.0; d = PetscRealPart(xit)+mu;
345: }
346: PetscCall(PEPGetST(pep,&st));
347: PetscCall(STGetMatStructure(st,&str));
348: PetscCall(TransformMatricesMoebius(A,str,a,b,c,d,T));
349: PetscFunctionReturn(PETSC_SUCCESS);
350: }
352: /* ------------------------------------------------------------------- */
353: /*
354: Auxiliary function to compute the residual norm of an eigenpair of a QEP defined
355: by coefficient matrices A
356: */
357: static PetscErrorCode PEPResidualNorm(Mat *A,PetscScalar kr,PetscScalar ki,Vec xr,Vec xi,Vec *z,PetscReal *norm)
358: {
359: PetscInt i,nmat=3;
360: PetscScalar vals[3];
361: Vec u,w;
362: #if !defined(PETSC_USE_COMPLEX)
363: Vec ui,wi;
364: PetscReal ni;
365: PetscBool imag;
366: PetscScalar ivals[3];
367: #endif
369: PetscFunctionBegin;
370: u = z[0]; w = z[1];
371: PetscCall(VecSet(u,0.0));
372: #if !defined(PETSC_USE_COMPLEX)
373: ui = z[2]; wi = z[3];
374: #endif
375: vals[0] = 1.0;
376: vals[1] = kr;
377: vals[2] = kr*kr-ki*ki;
378: #if !defined(PETSC_USE_COMPLEX)
379: ivals[0] = 0.0;
380: ivals[1] = ki;
381: ivals[2] = 2.0*kr*ki;
382: if (ki == 0 || PetscAbsScalar(ki) < PetscAbsScalar(kr*PETSC_MACHINE_EPSILON))
383: imag = PETSC_FALSE;
384: else {
385: imag = PETSC_TRUE;
386: PetscCall(VecSet(ui,0.0));
387: }
388: #endif
389: for (i=0;i<nmat;i++) {
390: if (vals[i]!=0.0) {
391: PetscCall(MatMult(A[i],xr,w));
392: PetscCall(VecAXPY(u,vals[i],w));
393: }
394: #if !defined(PETSC_USE_COMPLEX)
395: if (imag) {
396: if (ivals[i]!=0 || vals[i]!=0) {
397: PetscCall(MatMult(A[i],xi,wi));
398: if (vals[i]==0) PetscCall(MatMult(A[i],xr,w));
399: }
400: if (ivals[i]!=0) {
401: PetscCall(VecAXPY(u,-ivals[i],wi));
402: PetscCall(VecAXPY(ui,ivals[i],w));
403: }
404: if (vals[i]!=0) PetscCall(VecAXPY(ui,vals[i],wi));
405: }
406: #endif
407: }
408: PetscCall(VecNorm(u,NORM_2,norm));
409: #if !defined(PETSC_USE_COMPLEX)
410: if (imag) {
411: PetscCall(VecNorm(ui,NORM_2,&ni));
412: *norm = SlepcAbsEigenvalue(*norm,ni);
413: }
414: #endif
415: PetscFunctionReturn(PETSC_SUCCESS);
416: }
418: /* ------------------------------------------------------------------- */
419: /*
420: QEPDefiniteCheckError - check and print the residual norm of a transformed PEP
422: Input:
423: A: coefficient matrices of the original problem
424: pep: solver containing the computed solution of the transformed problem
425: xi,mu,hyperbolic: parameters used in transformation
426: */
427: PetscErrorCode QEPDefiniteCheckError(Mat *A,PEP pep,PetscBool hyperbolic,PetscReal xi,PetscReal mu)
428: {
429: PetscScalar er,ei;
430: PetscReal re,im,error;
431: Vec vr,vi,w[4];
432: PetscInt i,nconv;
433: BV bv;
434: char ex[30],sep[]=" ---------------------- --------------------\n";
436: PetscFunctionBegin;
437: PetscCall(PetscSNPrintf(ex,sizeof(ex),"||P(k)x||/||kx||"));
438: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%s k %s\n%s",sep,ex,sep));
439: PetscCall(PEPGetConverged(pep,&nconv));
440: PetscCall(PEPGetBV(pep,&bv));
441: PetscCall(BVCreateVec(bv,w));
442: PetscCall(VecDuplicate(w[0],&vr));
443: PetscCall(VecDuplicate(w[0],&vi));
444: for (i=1;i<4;i++) PetscCall(VecDuplicate(w[0],w+i));
445: for (i=0;i<nconv;i++) {
446: PetscCall(PEPGetEigenpair(pep,i,&er,&ei,vr,vi));
447: PetscCall(QEPDefiniteTransformMap(hyperbolic,xi,mu,1,&er,PETSC_TRUE));
448: PetscCall(PEPResidualNorm(A,er,0.0,vr,vi,w,&error));
449: error /= SlepcAbsEigenvalue(er,0.0);
450: #if defined(PETSC_USE_COMPLEX)
451: re = PetscRealPart(er);
452: im = PetscImaginaryPart(ei);
453: #else
454: re = er;
455: im = ei;
456: #endif
457: if (im!=0.0) PetscCall(PetscPrintf(PETSC_COMM_WORLD," % 9f%+9fi %12g\n",(double)re,(double)im,(double)error));
458: else PetscCall(PetscPrintf(PETSC_COMM_WORLD," % 12f %12g\n",(double)re,(double)error));
459: }
460: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%s",sep));
461: for (i=0;i<4;i++) PetscCall(VecDestroy(w+i));
462: PetscCall(VecDestroy(&vi));
463: PetscCall(VecDestroy(&vr));
464: PetscFunctionReturn(PETSC_SUCCESS);
465: }
467: /*TEST
469: testset:
470: requires: !single
471: args: -pep_nev 3 -nonhyperbolic -pep_target 2
472: output_file: output/ex40_1.out
473: filter: grep -v "Definite" | sed -e "s/iterations [0-9]\([0-9]*\)/iterations xx/g" | sed -e "s/[0-9]\.[0-9]*e[+-]\([0-9]*\)/removed/g"
474: test:
475: suffix: 1
476: requires: !complex
477: test:
478: suffix: 1_complex
479: requires: complex !mumps
480: test:
481: suffix: 1_transform
482: requires: !complex
483: args: -transform
484: test:
485: suffix: 1_transform_complex
486: requires: complex !mumps
487: args: -transform
489: TEST*/