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:   PetscObjectSetName((PetscObject)dms, "Particles");
 23:   DMSwarmInitializeFieldRegister(dms);
 24:   DMSwarmRegisterPetscDatatypeField(dms, "viscosity", 1, PETSC_REAL);
 25:   DMSwarmRegisterPetscDatatypeField(dms, "strain", 3, PETSC_REAL);
 26:   DMSwarmFinalizeFieldRegister(dms);
 27:   DMSwarmSetLocalSizes(dms, 5 + rank, 4);
 28:   DMView(dms, PETSC_VIEWER_STDOUT_WORLD);
 29:   DMSwarmGetLocalSize(dms, &nlocal);

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

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

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

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

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

 65: int main(int argc, char **argv)
 66: {
 68:   PetscInitialize(&argc, &argv, (char *)0, help);
 69:   ex2_1();
 70:   PetscFinalize();
 71:   return 0;
 72: }

 74: /*TEST

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

 82: TEST*/