Actual source code: ex5.c

  2: static char help[] = "Tests AODataRemap(). \n\n";

 4:  #include petscao.h

  8: int main(int argc,char **argv)
  9: {
 10:   int         n,nglobal,bs = 1,*keys,*data,ierr,rank,size,i,start,*news;
 11:   AOData      aodata;
 12:   AO          ao;

 14:   PetscInitialize(&argc,&argv,(char*)0,help);
 15:   PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);

 17:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank); n = rank + 2;
 18:   MPI_Allreduce(&n,&nglobal,1,MPI_INT,MPI_SUM,PETSC_COMM_WORLD);
 19:   MPI_Comm_size(PETSC_COMM_WORLD,&size);

 21:   /*
 22:        Create a database with one  key and one segment
 23:   */
 24:   AODataCreateBasic(PETSC_COMM_WORLD,&aodata);

 26:   /*
 27:        Put one segment in the key
 28:   */
 29:   AODataKeyAdd(aodata,"key1",PETSC_DECIDE,nglobal);

 31:   /* allocate space for the keys each processor will provide */
 32:   PetscMalloc(n*sizeof(int),&keys);

 34:   /*
 35:      We assign the first set of keys (0 to 2) to processor 0, etc.
 36:      This computes the first local key on each processor
 37:   */
 38:   MPI_Scan(&n,&start,1,MPI_INT,MPI_SUM,PETSC_COMM_WORLD);
 39:   start -= n;

 41:   for (i=0; i<n; i++) {
 42:     keys[i]     = start + i;
 43:   }

 45:   /* 
 46:       Allocate data for the first key and first segment 
 47:   */
 48:   PetscMalloc(bs*n*sizeof(int),&data);
 49:   for (i=0; i<n; i++) {
 50:     data[i]   = start + i + 1; /* the data is the neighbor to the right */
 51:   }
 52:   data[n-1] = 0; /* make it periodic */
 53:   AODataSegmentAdd(aodata,"key1","key1",bs,n,keys,data,PETSC_INT);
 54:   PetscFree(data);
 55:   PetscFree(keys);

 57:   /*
 58:         View the database
 59:   */
 60:   AODataView(aodata,PETSC_VIEWER_STDOUT_WORLD);
 61: 
 62:   /*
 63:          Remap the database so that i -> nglobal - i - 1
 64:   */
 65:   PetscMalloc(n*sizeof(int),&news);
 66:   for (i=0; i<n; i++) {
 67:     news[i] = nglobal - i - start - 1;
 68:   }
 69:   AOCreateBasic(PETSC_COMM_WORLD,n,news,PETSC_NULL,&ao);
 70:   PetscFree(news);
 71:   AODataKeyRemap(aodata,"key1",ao);
 72:   AODestroy(ao);
 73:   AODataView(aodata,PETSC_VIEWER_STDOUT_WORLD);
 74:   AODataDestroy(aodata);

 76:   PetscFinalize();
 77:   return 0;
 78: }
 79: