Actual source code: ex5.c

  1: /*$Id: ex5.c,v 1.19 2001/04/10 19:37:18 bsmith Exp $*/

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

 5:  #include petscao.h

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

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

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

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

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

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

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

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

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

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

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