1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests PetscOptionsGetScalar(), PetscOptionsScalarArray() for complex numbers\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscsys.h> 5c4762a1bSJed Brown 6c4762a1bSJed Brown int main(int argc,char **argv) 7c4762a1bSJed Brown { 8d0609cedSBarry Smith PetscInt n,i; 9c4762a1bSJed Brown PetscScalar a,array[10]; 10c4762a1bSJed Brown PetscReal rarray[10]; 11c4762a1bSJed Brown 12*327415f7SBarry Smith PetscFunctionBeginUser; 139566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); 149566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL,NULL,"-a",&a,NULL)); 159566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Scalar a = %g + %gi\n",(double)PetscRealPart(a),(double)PetscImaginaryPart(a))); 16c4762a1bSJed Brown 17d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"test options",NULL); 18c4762a1bSJed Brown n = 10; /* max num of input values */ 199566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-rarray", "Input a real array", "ex14.c", rarray, &n, NULL)); 20c4762a1bSJed Brown if (n) { 219566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Real rarray of length %" PetscInt_FMT "\n",n)); 22c4762a1bSJed Brown for (i=0; i<n; i++) { 239566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF," %g,\n",(double)rarray[i])); 24c4762a1bSJed Brown } 25c4762a1bSJed Brown } 26c4762a1bSJed Brown 27c4762a1bSJed Brown n = 10; /* max num of input values */ 289566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalarArray("-array", "Input a scalar array", "ex14.c", array, &n, NULL)); 29c4762a1bSJed Brown if (n) { 309566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Scalar rarray of length %" PetscInt_FMT "\n",n)); 31c4762a1bSJed Brown for (i=0; i<n; i++) { 32c4762a1bSJed Brown if (PetscImaginaryPart(array[i]) < 0.0) { 339566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF," %g - %gi\n",(double)PetscRealPart(array[i]),(double)PetscAbsReal(PetscImaginaryPart(array[i])))); 34c4762a1bSJed Brown } else { 359566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF," %g + %gi\n",(double)PetscRealPart(array[i]),(double)PetscImaginaryPart(array[i]))); 36c4762a1bSJed Brown } 37c4762a1bSJed Brown } 38c4762a1bSJed Brown } 39d0609cedSBarry Smith PetscOptionsEnd(); 409566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 41b122ec5aSJacob Faibussowitsch return 0; 42c4762a1bSJed Brown } 43c4762a1bSJed Brown 44c4762a1bSJed Brown /*TEST 45c4762a1bSJed Brown 46c4762a1bSJed Brown test: 47c4762a1bSJed Brown requires: complex 48c4762a1bSJed 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 49c4762a1bSJed Brown 50c4762a1bSJed Brown TEST*/ 51