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;
 14:   Vec            x;
 15:   PetscMPIInt    rank;
 16:   PetscInt       p,bs,nlocal;

 18:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);

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

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

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

 50:   DMSwarmCreateGlobalVectorFromField(dms,"viscosity",&x);
 51:   VecView(x,PETSC_VIEWER_STDOUT_WORLD);
 52:   DMSwarmDestroyGlobalVectorFromField(dms,"viscosity",&x);

 54:   DMSwarmCreateGlobalVectorFromField(dms,"strain",&x);
 55:   VecView(x,PETSC_VIEWER_STDOUT_WORLD);
 56:   DMSwarmDestroyGlobalVectorFromField(dms,"strain",&x);

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

 66: int main(int argc,char **argv)
 67: {

 69:   PetscInitialize(&argc,&argv,(char*)0,help);
 70:   ex2_1();
 71:   PetscFinalize();
 72:   return 0;
 73: }

 75: /*TEST

 77:    test:
 78:       requires: !complex double
 79:       nsize: 3
 80:       filter: grep -v atomic
 81:       filter_output: grep -v atomic

 83: TEST*/