Actual source code: ex1.c

petsc-3.7.5 2017-01-01
Report Typos and Errors
  2: static char help[] = "Tests PetscRandom functions.\n\n";

  4: #include <petscsys.h>

  6: /* Usage:
  7:    mpiexec -n <np> ./ex1 -n <num_of_random_numbers> -random_type <type> -log_summary
  8:                          -view_randomvalues <view_rank>
  9:                          -random_view ascii -random_view :filename
 10: */

 14: int main(int argc,char **argv)
 15: {
 16:   PetscInt       i,n = 1000,*values;
 17:   PetscRandom    rnd;
 18:   PetscScalar    value,avg = 0.0;
 20:   PetscMPIInt    rank;
 21:   PetscInt       view_rank=-1;
 22: #if defined(PETSC_USE_LOG)
 23:   PetscLogEvent  event;
 24: #endif

 26:   PetscInitialize(&argc,&argv,(char*)0,help);
 27:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
 28:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 29:   PetscOptionsGetInt(NULL,NULL,"-view_randomvalues",&view_rank,NULL);

 31:   PetscRandomCreate(PETSC_COMM_WORLD,&rnd);
 32:   /* force imaginary part of random number to always be zero; thus obtain reproducible results with real and complex numbers */
 33:   PetscRandomSetInterval(rnd,0.0,1.0);
 34:   PetscRandomSetFromOptions(rnd);

 36:   PetscMalloc1(n,&values);
 37:   for (i=0; i<n; i++) {
 38:     PetscRandomGetValue(rnd,&value);
 39:     avg += value;
 40:     if (view_rank == (PetscInt)rank) {
 41:       PetscPrintf(PETSC_COMM_SELF,"[%d] value[%D] = %6.4e\n",rank,i,(double)PetscRealPart(value));
 42:     }
 43:     values[i] = (PetscInt)(n*PetscRealPart(value) + 2.0);
 44:   }
 45:   avg = avg/((PetscReal)n);
 46:   if (view_rank == (PetscInt)rank) {
 47:     PetscPrintf(PETSC_COMM_SELF,"[%d] Average value %6.4e\n",rank,(double)PetscRealPart(avg));
 48:   }

 50:   PetscSortInt(n,values);

 52:   PetscLogEventRegister("Sort",0,&event);
 53:   PetscLogEventBegin(event,0,0,0,0);

 55:   PetscRandomSeed(rnd);
 56:   for (i=0; i<n; i++) {
 57:     PetscRandomGetValue(rnd,&value);
 58:     values[i] = (PetscInt)(n*PetscRealPart(value) + 2.0);
 59:     /* printf("value[%d] = %g\n",i,value); */
 60:   }
 61:   PetscSortInt(n,values);
 62:   PetscLogEventEnd(event,0,0,0,0);

 64:   for (i=1; i<n; i++) {
 65:     if (values[i] < values[i-1]) SETERRQ(PETSC_COMM_SELF,1,"Values not sorted");
 66:   }
 67:   PetscFree(values);
 68:   PetscRandomDestroy(&rnd);

 70:   PetscFinalize();
 71:   return 0;
 72: }