Actual source code: slda.c

  1: #include <../src/ts/characteristic/impls/da/slda.h>       /*I  "petsccharacteristic.h"  I*/

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

 12:   /* Pull out field names from DM */
 13:   PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);
 14:   PetscTypeCompare((PetscObject) viewer, PETSCVIEWERSTRING, &isstring);
 15:   if (iascii) {
 16:     PetscViewerASCIIPrintf(viewer,"  DMDA: dummy=%D\n", da->dummy);
 17:   } else if (isstring) {
 18:     PetscViewerStringSPrintf(viewer,"dummy %D", da->dummy);
 19:   } else {
 20:     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP, "Viewer type %s not supported for Characteristic DMDA", ((PetscObject) viewer)->type_name);
 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_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:   PetscMalloc(c->queueMax * sizeof(CharacteristicPointDA2D), &c->queue);
 62:   c->queueSize = 0;

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

 80: PetscErrorCode CharacteristicCreate_DA(Characteristic c)
 81: {
 82:   Characteristic_DA *da;
 83:   PetscErrorCode     ierr;

 86:   PetscNew(Characteristic_DA, &da);
 87:   PetscMemzero(da, sizeof(Characteristic_DA));
 88:   PetscLogObjectMemory(c, sizeof(Characteristic_DA));
 89:   c->data = (void *) da;

 91:   c->ops->setup   = CharacteristicSetUp_DA;
 92:   c->ops->destroy = CharacteristicDestroy_DA;
 93:   c->ops->view    = CharacteristicView_DA;

 95:   da->dummy = 0;
 96:   return(0);
 97: }

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

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

115:   if (bx == DMDA_BOUNDARY_PERIODIC) {
116:       while (*x >= ( PetscScalar ) gx ) { *x -= ( PetscScalar ) gx; }
117:       while (*x < 0.0 )                 { *x += ( PetscScalar ) gx; }
118:     }
119:     if (by == DMDA_BOUNDARY_PERIODIC) {
120:       while (*y >= ( PetscScalar ) gy ) { *y -= ( PetscScalar ) gy; }
121:       while (*y < 0.0 )                 { *y += ( PetscScalar ) gy; }
122:     }
123:   PetscFunctionReturn(ierr);
124: }