xref: /petsc/src/sys/tests/ex14.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown 
2*c4762a1bSJed Brown static char help[] = "Tests PetscOptionsGetScalar(), PetscOptionsScalarArray() for complex numbers\n";
3*c4762a1bSJed Brown 
4*c4762a1bSJed Brown #include <petscsys.h>
5*c4762a1bSJed Brown 
6*c4762a1bSJed Brown int main(int argc,char **argv)
7*c4762a1bSJed Brown {
8*c4762a1bSJed Brown   PetscInt    ierr,n,i;
9*c4762a1bSJed Brown   PetscScalar a,array[10];
10*c4762a1bSJed Brown   PetscReal   rarray[10];
11*c4762a1bSJed Brown 
12*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
13*c4762a1bSJed Brown   ierr = PetscOptionsGetScalar(NULL,NULL,"-a",&a,NULL);CHKERRQ(ierr);
14*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_SELF,"Scalar a = %g + %gi\n",(double)PetscRealPart(a),(double)PetscImaginaryPart(a));CHKERRQ(ierr);
15*c4762a1bSJed Brown 
16*c4762a1bSJed Brown   ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"test options",NULL);CHKERRQ(ierr);
17*c4762a1bSJed Brown   n = 10; /* max num of input values */
18*c4762a1bSJed Brown   ierr = PetscOptionsRealArray("-rarray", "Input a real array", "ex14.c", rarray, &n, NULL);CHKERRQ(ierr);
19*c4762a1bSJed Brown   if (n) {
20*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_SELF,"Real rarray of length %d\n",n);CHKERRQ(ierr);
21*c4762a1bSJed Brown     for (i=0; i<n; i++){
22*c4762a1bSJed Brown       ierr = PetscPrintf(PETSC_COMM_SELF," %g,\n",rarray[i]);CHKERRQ(ierr);
23*c4762a1bSJed Brown     }
24*c4762a1bSJed Brown   }
25*c4762a1bSJed Brown 
26*c4762a1bSJed Brown   n = 10; /* max num of input values */
27*c4762a1bSJed Brown   ierr = PetscOptionsScalarArray("-array", "Input a scalar array", "ex14.c", array, &n, NULL);CHKERRQ(ierr);
28*c4762a1bSJed Brown   if (n) {
29*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_SELF,"Scalar rarray of length %d\n",n);CHKERRQ(ierr);
30*c4762a1bSJed Brown     for (i=0; i<n; i++){
31*c4762a1bSJed Brown       if (PetscImaginaryPart(array[i]) < 0.0) {
32*c4762a1bSJed Brown         ierr = PetscPrintf(PETSC_COMM_SELF," %g - %gi\n",(double)PetscRealPart(array[i]),(double)PetscAbsReal(PetscImaginaryPart(array[i])));CHKERRQ(ierr);
33*c4762a1bSJed Brown       } else {
34*c4762a1bSJed Brown         ierr = PetscPrintf(PETSC_COMM_SELF," %g + %gi\n",(double)PetscRealPart(array[i]),(double)PetscImaginaryPart(array[i]));CHKERRQ(ierr);
35*c4762a1bSJed Brown       }
36*c4762a1bSJed Brown     }
37*c4762a1bSJed Brown   }
38*c4762a1bSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
39*c4762a1bSJed Brown   ierr = PetscFinalize();
40*c4762a1bSJed Brown   return ierr;
41*c4762a1bSJed Brown }
42*c4762a1bSJed Brown 
43*c4762a1bSJed Brown 
44*c4762a1bSJed Brown 
45*c4762a1bSJed Brown /*TEST
46*c4762a1bSJed Brown 
47*c4762a1bSJed Brown    test:
48*c4762a1bSJed Brown       requires: complex
49*c4762a1bSJed Brown       args: -array 1.0,-2-3i,4.5+6.2i,4.5,6.8+4i,i,-i,-1.2i -rarray 1,2,3 -a 1.5+2.1i
50*c4762a1bSJed Brown 
51*c4762a1bSJed Brown TEST*/
52