Actual source code: slda.c
1: #include "slda.h" /*I "characteristic.h" I*/
5: PetscErrorCode CharacteristicView_DA(Characteristic c, PetscViewer viewer)
6: {
7: Characteristic_DA *da = (Characteristic_DA *) c->data;
8: PetscTruth iascii, isstring;
9: PetscErrorCode ierr;
12: /* Pull out field names from DA */
13: PetscTypeCompare((PetscObject) viewer, PETSC_VIEWER_ASCII, &iascii);
14: PetscTypeCompare((PetscObject) viewer, PETSC_VIEWER_STRING, &isstring);
15: if (iascii) {
16: PetscViewerASCIIPrintf(viewer," DA: dummy=%D\n", da->dummy);
17: } else if (isstring) {
18: PetscViewerStringSPrintf(viewer,"dummy %D", da->dummy);
19: } else {
20: SETERRQ1(PETSC_ERR_SUP, "Viewer type %s not supported for Characteristic DA", ((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: DAGetInfo(c->velocityDA, &dim, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0);
48: if (c->structured) {
49: c->numIds = dim;
50: } else {
51: c->numIds = 3;
52: }
53: if (c->numFieldComp > MAX_COMPONENTS) {
54: SETERRQ2(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);
55: }
56: numValues = 4 + MAX_COMPONENTS;
58: /* Create new MPI datatype for communication of characteristic point structs */
59: blockLen[0] = 1+c->numIds; indices[0] = 0; oldtypes[0] = MPIU_INT;
60: blockLen[1] = numValues; indices[1] = (1+c->numIds)*sizeof(PetscInt); oldtypes[1] = MPIU_SCALAR;
61: MPI_Type_struct(2, blockLen, indices, oldtypes, &c->itemType);
62: MPI_Type_commit(&c->itemType);
64: /* Initialize the local queue for char foot values */
65: VecGetLocalSize(c->velocity, &c->queueMax);
66: PetscMalloc(c->queueMax * sizeof(CharacteristicPointDA2D), &c->queue);
67: c->queueSize = 0;
69: /* Allocate communication structures */
70: if (c->numNeighbors <= 0) {
71: SETERRQ1(PETSC_ERR_ARG_WRONGSTATE, "Invalid number of neighbors %d. Call CharactersiticSetNeighbors() before setup.", c->numNeighbors);
72: }
73: PetscMalloc(c->numNeighbors * sizeof(PetscInt), &c->needCount);
74: PetscMalloc(c->numNeighbors * sizeof(PetscInt), &c->localOffsets);
75: PetscMalloc(c->numNeighbors * sizeof(PetscInt), &c->fillCount);
76: PetscMalloc(c->numNeighbors * sizeof(PetscInt), &c->remoteOffsets);
77: PetscMalloc((c->numNeighbors-1) * sizeof(MPI_Request), &c->request);
78: PetscMalloc((c->numNeighbors-1) * sizeof(MPI_Status), &c->status);
79: return(0);
80: }
85: PetscErrorCode CharacteristicCreate_DA(Characteristic c)
86: {
87: Characteristic_DA *da;
88: PetscErrorCode ierr;
91: PetscNew(Characteristic_DA, &da);
92: PetscMemzero(da, sizeof(Characteristic_DA));
93: PetscLogObjectMemory(c, sizeof(Characteristic_DA));
94: c->data = (void *) da;
96: c->ops->setup = CharacteristicSetUp_DA;
97: c->ops->destroy = CharacteristicDestroy_DA;
98: c->ops->view = CharacteristicView_DA;
100: da->dummy = 0;
101: return(0);
102: }
107: /* -----------------------------------------------------------------------------
108: Checks for periodicity of a DA and Maps points outside of a domain back onto the domain
109: using appropriate periodicity. At the moment assumes only a 2-D DA.
110: ----------------------------------------------------------------------------------------*/
111: PetscErrorCode DAMapCoordsToPeriodicDomain(DA da, PetscScalar *x, PetscScalar *y)
112: {
113: DAPeriodicType periodic_type;
114: PetscInt dim, gx, gy;
118: DAGetInfo(da, &dim, &gx, &gy, 0, 0, 0, 0, 0, 0, &periodic_type, 0);
120: if ( periodic_type == DA_NONPERIODIC ) {
121: 0;
122: } else {
123: if (periodic_type==DA_XPERIODIC || periodic_type==DA_XYPERIODIC) {
124: while (*x >= ( PetscScalar ) gx ) { *x -= ( PetscScalar ) gx; }
125: while (*x < 0.0 ) { *x += ( PetscScalar ) gx; }
126: }
127: if (periodic_type==DA_YPERIODIC || periodic_type==DA_XYPERIODIC) {
128: while (*y >= ( PetscScalar ) gy ) { *y -= ( PetscScalar ) gy; }
129: while (*y < 0.0 ) { *y += ( PetscScalar ) gy; }
130: }
131: }
132:
133: PetscFunctionReturn(ierr);
134: }