Actual source code: swarm_ex2.c

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

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

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

 17:   PetscFunctionBegin;
 18:   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));

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

 31:   {
 32:     PetscReal *array;
 33:     PetscCall(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:     PetscCall(DMSwarmRestoreField(dms, "viscosity", &bs, NULL, (void **)&array));
 36:   }

 38:   {
 39:     PetscReal *array;
 40:     PetscCall(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:     PetscCall(DMSwarmRestoreField(dms, "strain", &bs, NULL, (void **)&array));
 47:   }

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

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

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

 65: int main(int argc, char **argv)
 66: {
 67:   PetscFunctionBeginUser;
 68:   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
 69:   PetscCall(ex2_1());
 70:   PetscCall(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*/