Actual source code: slda.c

petsc-3.5.4 2015-05-23
Report Typos and Errors
  1: #include <../src/ts/characteristic/impls/da/slda.h>       /*I  "petsccharacteristic.h"  I*/
  2: #include <petscdmda.h>
  3: #include <petscviewer.h>

  7: PetscErrorCode CharacteristicView_DA(Characteristic c, PetscViewer viewer)
  8: {
  9:   Characteristic_DA *da = (Characteristic_DA*) c->data;
 10:   PetscBool         iascii, isstring;
 11:   PetscErrorCode    ierr;

 14:   /* Pull out field names from DM */
 15:   PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);
 16:   PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERSTRING, &isstring);
 17:   if (iascii) {
 18:     PetscViewerASCIIPrintf(viewer,"  DMDA: dummy=%D\n", da->dummy);
 19:   } else if (isstring) {
 20:     PetscViewerStringSPrintf(viewer,"dummy %D", da->dummy);
 21:   }
 22:   return(0);
 23: }

 27: PetscErrorCode CharacteristicDestroy_DA(Characteristic c)
 28: {
 29:   Characteristic_DA *da = (Characteristic_DA*) c->data;
 30:   PetscErrorCode    ierr;

 33:   PetscFree(da);
 34:   return(0);
 35: }

 39: PetscErrorCode CharacteristicSetUp_DA(Characteristic c)
 40: {
 41:   PetscMPIInt    blockLen[2];
 42:   MPI_Aint       indices[2];
 43:   MPI_Datatype   oldtypes[2];
 44:   PetscInt       dim, numValues;

 47:   DMDAGetInfo(c->velocityDA, &dim, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0,0,0);
 48:   if (c->structured) c->numIds = dim;
 49:   else c->numIds = 3;
 50:   if (c->numFieldComp > MAX_COMPONENTS) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE, "The maximum number of fields allowed is %d, you have %d. You can recompile after increasing MAX_COMPONENTS.", MAX_COMPONENTS, c->numFieldComp);
 51:   numValues = 4 + MAX_COMPONENTS;

 53:   /* Create new MPI datatype for communication of characteristic point structs */
 54:   blockLen[0] = 1+c->numIds; indices[0] = 0;                              oldtypes[0] = MPIU_INT;
 55:   blockLen[1] = numValues;   indices[1] = (1+c->numIds)*sizeof(PetscInt); oldtypes[1] = MPIU_SCALAR;
 56:   MPI_Type_create_struct(2, blockLen, indices, oldtypes, &c->itemType);
 57:   MPI_Type_commit(&c->itemType);

 59:   /* Initialize the local queue for char foot values */
 60:   VecGetLocalSize(c->velocity, &c->queueMax);
 61:   PetscMalloc1(c->queueMax, &c->queue);
 62:   c->queueSize = 0;

 64:   /* Allocate communication structures */
 65:   if (c->numNeighbors <= 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Invalid number of neighbors %d. Call CharactersiticSetNeighbors() before setup.", c->numNeighbors);
 66:   PetscMalloc1(c->numNeighbors, &c->needCount);
 67:   PetscMalloc1(c->numNeighbors, &c->localOffsets);
 68:   PetscMalloc1(c->numNeighbors, &c->fillCount);
 69:   PetscMalloc1(c->numNeighbors, &c->remoteOffsets);
 70:   PetscMalloc1((c->numNeighbors-1), &c->request);
 71:   PetscMalloc1((c->numNeighbors-1),  &c->status);
 72:   return(0);
 73: }

 77: PETSC_EXTERN PetscErrorCode CharacteristicCreate_DA(Characteristic c)
 78: {
 79:   Characteristic_DA *da;
 80:   PetscErrorCode    ierr;

 83:   PetscNew(&da);
 84:   PetscMemzero(da, sizeof(Characteristic_DA));
 85:   PetscLogObjectMemory((PetscObject)c, sizeof(Characteristic_DA));
 86:   c->data = (void*) da;

 88:   c->ops->setup   = CharacteristicSetUp_DA;
 89:   c->ops->destroy = CharacteristicDestroy_DA;
 90:   c->ops->view    = CharacteristicView_DA;

 92:   da->dummy = 0;
 93:   return(0);
 94: }

 98: /* -----------------------------------------------------------------------------
 99:    Checks for periodicity of a DM and Maps points outside of a domain back onto the domain
100:    using appropriate periodicity. At the moment assumes only a 2-D DMDA.
101:    ----------------------------------------------------------------------------------------*/
102: PetscErrorCode DMDAMapCoordsToPeriodicDomain(DM da, PetscScalar *x, PetscScalar *y)
103: {
104:   DMBoundaryType bx, by;
105:   PetscInt       dim, gx, gy;

109:   DMDAGetInfo(da, &dim, &gx, &gy, 0, 0, 0, 0, 0, 0, &bx, &by, 0, 0);

111:   if (bx == DM_BOUNDARY_PERIODIC) {
112:       while (*x >= (PetscScalar)gx) *x -= (PetscScalar)gx;
113:       while (*x < 0.0)              *x += (PetscScalar)gx;
114:     }
115:     if (by == DM_BOUNDARY_PERIODIC) {
116:       while (*y >= (PetscScalar)gy) *y -= (PetscScalar)gy;
117:       while (*y < 0.0)              *y += (PetscScalar)gy;
118:     }
119:   PetscFunctionReturn(ierr);
120: }