Actual source code: ex6.c

  1: /*$Id: ex6.c,v 1.20 2001/08/10 03:34:29 bsmith Exp $*/

  3: static char help[] = "Tests removing entries from an AOData. \n\n";

 5:  #include petscao.h

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

 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 two sets of keys 
 24:   */
 25:   AODataCreateBasic(PETSC_COMM_WORLD,&aodata);

 27:   /*
 28:        Put two segments in the first key and one in the second
 29:   */
 30:   AODataKeyAdd(aodata,"key1",PETSC_DECIDE,nglobal);
 31:   AODataKeyAdd(aodata,"key2",PETSC_DECIDE,nglobal);

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

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

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

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

 58:   /*
 59:       Allocate data for first key and second segment 
 60:   */
 61:   bs   = 3;
 62:   PetscMalloc(bs*n*sizeof(PetscReal),&gd);
 63:   for (i=0; i<n; i++) {
 64:     gd[3*i]   = -(start + i);
 65:     gd[3*i+1] = -(start + i) - 10000;
 66:     gd[3*i+2] = -(start + i) - 100000;
 67:   }
 68:   AODataSegmentAdd(aodata,"key1","seg2",bs,n,keys,gd,PETSC_REAL);

 70:   /*
 71:        Use same data for second key and first segment 
 72:   */
 73:   AODataSegmentAdd(aodata,"key2","seg1",bs,n,keys,gd,PETSC_REAL);
 74:   PetscFree(gd);
 75:   PetscFree(keys);

 77:   /*
 78:      View the database
 79:   */
 80:   AODataView(aodata,PETSC_VIEWER_STDOUT_WORLD);

 82:   /*
 83:        Remove a key and a single segment from the database
 84:   */
 85:   AODataKeyRemove(aodata,"key2");
 86:   AODataSegmentRemove(aodata,"key1","seg1");

 88:   AODataView(aodata,PETSC_VIEWER_STDOUT_WORLD);

 90:   AODataDestroy(aodata);

 92:   PetscFinalize();
 93:   return 0;
 94: }
 95: