Actual source code: swarm_ex2.c


  2: static char help[] = "Tests DMSwarm\n\n";

  4: #include <petscdm.h>
  5: #include <petscdmda.h>
  6: #include <petscdmswarm.h>

  8: /*
  9:  Checks for variable blocksize
 10: */
 11: PetscErrorCode ex2_1(void)
 12: {
 13:   DM             dms;
 15:   Vec            x;
 16:   PetscMPIInt    rank;
 17:   PetscInt       p,bs,nlocal;

 20:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);

 22:   DMCreate(PETSC_COMM_WORLD,&dms);
 23:   DMSetType(dms,DMSWARM);
 24:   DMSwarmInitializeFieldRegister(dms);
 25:   DMSwarmRegisterPetscDatatypeField(dms,"viscosity",1,PETSC_REAL);
 26:   DMSwarmRegisterPetscDatatypeField(dms,"strain",3,PETSC_REAL);
 27:   DMSwarmFinalizeFieldRegister(dms);
 28:   DMSwarmSetLocalSizes(dms,5+rank,4);
 29:   DMView(dms,PETSC_VIEWER_STDOUT_WORLD);
 30:   DMSwarmGetLocalSize(dms,&nlocal);

 32:   {
 33:     PetscReal *array;
 34:     DMSwarmGetField(dms,"viscosity",&bs,NULL,(void**)&array);
 35:     for (p=0; p<nlocal; p++) {
 36:       array[p] = 11.1 + p*0.1 + rank*100.0;
 37:     }
 38:     DMSwarmRestoreField(dms,"viscosity",&bs,NULL,(void**)&array);
 39:   }

 41:   {
 42:     PetscReal *array;
 43:     DMSwarmGetField(dms,"strain",&bs,NULL,(void**)&array);
 44:     for (p=0; p<nlocal; p++) {
 45:       array[bs*p+0] = 2.0e-2 + p*0.001 + rank*1.0;
 46:       array[bs*p+1] = 2.0e-2 + p*0.002 + rank*1.0;
 47:       array[bs*p+2] = 2.0e-2 + p*0.003 + rank*1.0;
 48:     }
 49:     DMSwarmRestoreField(dms,"strain",&bs,NULL,(void**)&array);
 50:   }

 52:   DMSwarmCreateGlobalVectorFromField(dms,"viscosity",&x);
 53:   VecView(x,PETSC_VIEWER_STDOUT_WORLD);
 54:   DMSwarmDestroyGlobalVectorFromField(dms,"viscosity",&x);

 56:   DMSwarmCreateGlobalVectorFromField(dms,"strain",&x);
 57:   VecView(x,PETSC_VIEWER_STDOUT_WORLD);
 58:   DMSwarmDestroyGlobalVectorFromField(dms,"strain",&x);

 60:   DMSwarmVectorDefineField(dms,"strain");
 61:   DMCreateGlobalVector(dms,&x);
 62:   VecView(x,PETSC_VIEWER_STDOUT_WORLD);
 63:   VecDestroy(&x);
 64:   DMDestroy(&dms);
 65:   return(0);
 66: }

 68: int main(int argc,char **argv)
 69: {

 72:   PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
 73:   ex2_1();
 74:   PetscFinalize();
 75:   return ierr;
 76: }

 78: /*TEST

 80:    test:
 81:       requires: !complex double
 82:       nsize: 3
 83:       filter: grep -v atomic
 84:       filter_output: grep -v atomic

 86: TEST*/