10e2ec84fSDave May #include <petscdm.h>
2cc4c1da9SBarry Smith #include <petscdmda.h> /*I "petscdmda.h" I*/
3cc4c1da9SBarry Smith #include <petsc/private/dmswarmimpl.h> /*I "petscdmswarm.h" I*/
4279f676cSBarry Smith #include "../src/dm/impls/swarm/data_bucket.h"
50e2ec84fSDave May
private_DMSwarmCreateCellLocalCoords_DA_Q1_Regular(PetscInt dim,PetscInt np[],PetscInt * _npoints,PetscReal ** _xi)666976f2fSJacob Faibussowitsch static PetscErrorCode private_DMSwarmCreateCellLocalCoords_DA_Q1_Regular(PetscInt dim, PetscInt np[], PetscInt *_npoints, PetscReal **_xi)
7d71ae5a4SJacob Faibussowitsch {
80e2ec84fSDave May PetscReal *xi;
9f9368e01SSatish Balay PetscInt d, npoints = 0, cnt;
100e2ec84fSDave May PetscReal ds[] = {0.0, 0.0, 0.0};
110e2ec84fSDave May PetscInt ii, jj, kk;
120e2ec84fSDave May
130e2ec84fSDave May PetscFunctionBegin;
140e2ec84fSDave May switch (dim) {
15d71ae5a4SJacob Faibussowitsch case 1:
16d71ae5a4SJacob Faibussowitsch npoints = np[0];
17d71ae5a4SJacob Faibussowitsch break;
18d71ae5a4SJacob Faibussowitsch case 2:
19d71ae5a4SJacob Faibussowitsch npoints = np[0] * np[1];
20d71ae5a4SJacob Faibussowitsch break;
21d71ae5a4SJacob Faibussowitsch case 3:
22d71ae5a4SJacob Faibussowitsch npoints = np[0] * np[1] * np[2];
23d71ae5a4SJacob Faibussowitsch break;
240e2ec84fSDave May }
25ad540459SPierre Jolivet for (d = 0; d < dim; d++) ds[d] = 2.0 / ((PetscReal)np[d]);
260e2ec84fSDave May
279566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dim * npoints, &xi));
280e2ec84fSDave May switch (dim) {
290e2ec84fSDave May case 1:
300e2ec84fSDave May cnt = 0;
310e2ec84fSDave May for (ii = 0; ii < np[0]; ii++) {
320e2ec84fSDave May xi[dim * cnt + 0] = -1.0 + 0.5 * ds[d] + ii * ds[0];
330e2ec84fSDave May cnt++;
340e2ec84fSDave May }
350e2ec84fSDave May break;
360e2ec84fSDave May
370e2ec84fSDave May case 2:
380e2ec84fSDave May cnt = 0;
390e2ec84fSDave May for (jj = 0; jj < np[1]; jj++) {
400e2ec84fSDave May for (ii = 0; ii < np[0]; ii++) {
410e2ec84fSDave May xi[dim * cnt + 0] = -1.0 + 0.5 * ds[0] + ii * ds[0];
420e2ec84fSDave May xi[dim * cnt + 1] = -1.0 + 0.5 * ds[1] + jj * ds[1];
430e2ec84fSDave May cnt++;
440e2ec84fSDave May }
450e2ec84fSDave May }
460e2ec84fSDave May break;
470e2ec84fSDave May
480e2ec84fSDave May case 3:
490e2ec84fSDave May cnt = 0;
500e2ec84fSDave May for (kk = 0; kk < np[2]; kk++) {
510e2ec84fSDave May for (jj = 0; jj < np[1]; jj++) {
520e2ec84fSDave May for (ii = 0; ii < np[0]; ii++) {
530e2ec84fSDave May xi[dim * cnt + 0] = -1.0 + 0.5 * ds[0] + ii * ds[0];
540e2ec84fSDave May xi[dim * cnt + 1] = -1.0 + 0.5 * ds[1] + jj * ds[1];
550e2ec84fSDave May xi[dim * cnt + 2] = -1.0 + 0.5 * ds[2] + kk * ds[2];
560e2ec84fSDave May cnt++;
570e2ec84fSDave May }
580e2ec84fSDave May }
590e2ec84fSDave May }
600e2ec84fSDave May break;
610e2ec84fSDave May }
620e2ec84fSDave May *_npoints = npoints;
630e2ec84fSDave May *_xi = xi;
643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
650e2ec84fSDave May }
660e2ec84fSDave May
private_DMSwarmCreateCellLocalCoords_DA_Q1_Gauss(PetscInt dim,PetscInt np_1d,PetscInt * _npoints,PetscReal ** _xi)6766976f2fSJacob Faibussowitsch static PetscErrorCode private_DMSwarmCreateCellLocalCoords_DA_Q1_Gauss(PetscInt dim, PetscInt np_1d, PetscInt *_npoints, PetscReal **_xi)
68d71ae5a4SJacob Faibussowitsch {
690e2ec84fSDave May PetscQuadrature quadrature;
700e2ec84fSDave May const PetscReal *quadrature_xi;
710e2ec84fSDave May PetscReal *xi;
720e2ec84fSDave May PetscInt d, q, npoints_q;
730e2ec84fSDave May
740e2ec84fSDave May PetscFunctionBegin;
759566063dSJacob Faibussowitsch PetscCall(PetscDTGaussTensorQuadrature(dim, 1, np_1d, -1.0, 1.0, &quadrature));
769566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(quadrature, NULL, NULL, &npoints_q, &quadrature_xi, NULL));
779566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dim * npoints_q, &xi));
780e2ec84fSDave May for (q = 0; q < npoints_q; q++) {
79ad540459SPierre Jolivet for (d = 0; d < dim; d++) xi[dim * q + d] = quadrature_xi[dim * q + d];
800e2ec84fSDave May }
819566063dSJacob Faibussowitsch PetscCall(PetscQuadratureDestroy(&quadrature));
820e2ec84fSDave May *_npoints = npoints_q;
830e2ec84fSDave May *_xi = xi;
843ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
850e2ec84fSDave May }
860e2ec84fSDave May
private_DMSwarmInsertPointsUsingCellDM_DA_Q1(DM dm,DM dmc,PetscInt npoints,DMSwarmPICLayoutType layout)8766976f2fSJacob Faibussowitsch static PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_DA_Q1(DM dm, DM dmc, PetscInt npoints, DMSwarmPICLayoutType layout)
88d71ae5a4SJacob Faibussowitsch {
89*19307e5cSMatthew G. Knepley DMSwarmCellDM celldm;
900e2ec84fSDave May PetscInt dim, npoints_q;
910e2ec84fSDave May PetscInt nel, npe, e, q, k, d;
920e2ec84fSDave May const PetscInt *element_list;
930e2ec84fSDave May PetscReal **basis;
940e2ec84fSDave May PetscReal *xi;
950e2ec84fSDave May Vec coor;
960e2ec84fSDave May const PetscScalar *_coor;
970e2ec84fSDave May PetscReal *elcoor;
980e2ec84fSDave May PetscReal *swarm_coor;
990e2ec84fSDave May PetscInt *swarm_cellid;
100*19307e5cSMatthew G. Knepley PetscInt pcnt, Nfc;
101*19307e5cSMatthew G. Knepley const char **coordFields, *cellid;
1020e2ec84fSDave May
1030e2ec84fSDave May PetscFunctionBegin;
1049566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim));
1050e2ec84fSDave May switch (layout) {
1069371c9d4SSatish Balay case DMSWARMPIC_LAYOUT_REGULAR: {
1070e2ec84fSDave May PetscInt np_dir[3];
1080e2ec84fSDave May np_dir[0] = np_dir[1] = np_dir[2] = npoints;
1099566063dSJacob Faibussowitsch PetscCall(private_DMSwarmCreateCellLocalCoords_DA_Q1_Regular(dim, np_dir, &npoints_q, &xi));
1109371c9d4SSatish Balay } break;
111d71ae5a4SJacob Faibussowitsch case DMSWARMPIC_LAYOUT_GAUSS:
112d71ae5a4SJacob Faibussowitsch PetscCall(private_DMSwarmCreateCellLocalCoords_DA_Q1_Gauss(dim, npoints, &npoints_q, &xi));
113d71ae5a4SJacob Faibussowitsch break;
114d3393ee1SDave May
1159371c9d4SSatish Balay case DMSWARMPIC_LAYOUT_SUBDIVISION: {
116d3393ee1SDave May PetscInt s, nsub;
117d3393ee1SDave May PetscInt np_dir[3];
118d3393ee1SDave May nsub = npoints;
119d3393ee1SDave May np_dir[0] = 1;
120ad540459SPierre Jolivet for (s = 0; s < nsub; s++) np_dir[0] *= 2;
121d3393ee1SDave May np_dir[1] = np_dir[0];
122d3393ee1SDave May np_dir[2] = np_dir[0];
1239566063dSJacob Faibussowitsch PetscCall(private_DMSwarmCreateCellLocalCoords_DA_Q1_Regular(dim, np_dir, &npoints_q, &xi));
1249371c9d4SSatish Balay } break;
125d71ae5a4SJacob Faibussowitsch default:
126d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "A valid DMSwarmPIC layout must be provided");
1270e2ec84fSDave May }
1280e2ec84fSDave May
1299566063dSJacob Faibussowitsch PetscCall(DMDAGetElements(dmc, &nel, &npe, &element_list));
1309566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dim * npe, &elcoor));
1319566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(npoints_q, &basis));
1320e2ec84fSDave May for (q = 0; q < npoints_q; q++) {
1339566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(npe, &basis[q]));
1340e2ec84fSDave May
1350e2ec84fSDave May switch (dim) {
1360e2ec84fSDave May case 1:
1370e2ec84fSDave May basis[q][0] = 0.5 * (1.0 - xi[dim * q + 0]);
1380e2ec84fSDave May basis[q][1] = 0.5 * (1.0 + xi[dim * q + 0]);
1390e2ec84fSDave May break;
1400e2ec84fSDave May case 2:
1410e2ec84fSDave May basis[q][0] = 0.25 * (1.0 - xi[dim * q + 0]) * (1.0 - xi[dim * q + 1]);
1420e2ec84fSDave May basis[q][1] = 0.25 * (1.0 + xi[dim * q + 0]) * (1.0 - xi[dim * q + 1]);
1430e2ec84fSDave May basis[q][2] = 0.25 * (1.0 + xi[dim * q + 0]) * (1.0 + xi[dim * q + 1]);
1440e2ec84fSDave May basis[q][3] = 0.25 * (1.0 - xi[dim * q + 0]) * (1.0 + xi[dim * q + 1]);
1450e2ec84fSDave May break;
1460e2ec84fSDave May
1470e2ec84fSDave May case 3:
1480e2ec84fSDave May basis[q][0] = 0.125 * (1.0 - xi[dim * q + 0]) * (1.0 - xi[dim * q + 1]) * (1.0 - xi[dim * q + 2]);
1490e2ec84fSDave May basis[q][1] = 0.125 * (1.0 + xi[dim * q + 0]) * (1.0 - xi[dim * q + 1]) * (1.0 - xi[dim * q + 2]);
1500e2ec84fSDave May basis[q][2] = 0.125 * (1.0 + xi[dim * q + 0]) * (1.0 + xi[dim * q + 1]) * (1.0 - xi[dim * q + 2]);
1510e2ec84fSDave May basis[q][3] = 0.125 * (1.0 - xi[dim * q + 0]) * (1.0 + xi[dim * q + 1]) * (1.0 - xi[dim * q + 2]);
1520e2ec84fSDave May basis[q][4] = 0.125 * (1.0 - xi[dim * q + 0]) * (1.0 - xi[dim * q + 1]) * (1.0 + xi[dim * q + 2]);
1530e2ec84fSDave May basis[q][5] = 0.125 * (1.0 + xi[dim * q + 0]) * (1.0 - xi[dim * q + 1]) * (1.0 + xi[dim * q + 2]);
1540e2ec84fSDave May basis[q][6] = 0.125 * (1.0 + xi[dim * q + 0]) * (1.0 + xi[dim * q + 1]) * (1.0 + xi[dim * q + 2]);
1550e2ec84fSDave May basis[q][7] = 0.125 * (1.0 - xi[dim * q + 0]) * (1.0 + xi[dim * q + 1]) * (1.0 + xi[dim * q + 2]);
1560e2ec84fSDave May break;
1570e2ec84fSDave May }
1580e2ec84fSDave May }
1590e2ec84fSDave May
160*19307e5cSMatthew G. Knepley PetscCall(DMSwarmGetCellDMActive(dm, &celldm));
161*19307e5cSMatthew G. Knepley PetscCall(DMSwarmCellDMGetCoordinateFields(celldm, &Nfc, &coordFields));
162*19307e5cSMatthew G. Knepley PetscCheck(Nfc == 1, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "We only support a single coordinate field right now, not %" PetscInt_FMT, Nfc);
163*19307e5cSMatthew G. Knepley PetscCall(DMSwarmCellDMGetCellID(celldm, &cellid));
164*19307e5cSMatthew G. Knepley
1659566063dSJacob Faibussowitsch PetscCall(DMSwarmSetLocalSizes(dm, npoints_q * nel, -1));
166*19307e5cSMatthew G. Knepley PetscCall(DMSwarmGetField(dm, coordFields[0], NULL, NULL, (void **)&swarm_coor));
167*19307e5cSMatthew G. Knepley PetscCall(DMSwarmGetField(dm, cellid, NULL, NULL, (void **)&swarm_cellid));
1680e2ec84fSDave May
1699566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dmc, &coor));
1709566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(coor, &_coor));
1710e2ec84fSDave May pcnt = 0;
1720e2ec84fSDave May for (e = 0; e < nel; e++) {
1730e2ec84fSDave May const PetscInt *element = &element_list[npe * e];
1740e2ec84fSDave May
1750e2ec84fSDave May for (k = 0; k < npe; k++) {
176ad540459SPierre Jolivet for (d = 0; d < dim; d++) elcoor[dim * k + d] = PetscRealPart(_coor[dim * element[k] + d]);
1770e2ec84fSDave May }
1780e2ec84fSDave May
1790e2ec84fSDave May for (q = 0; q < npoints_q; q++) {
180ad540459SPierre Jolivet for (d = 0; d < dim; d++) swarm_coor[dim * pcnt + d] = 0.0;
1810e2ec84fSDave May for (k = 0; k < npe; k++) {
182ad540459SPierre Jolivet for (d = 0; d < dim; d++) swarm_coor[dim * pcnt + d] += basis[q][k] * elcoor[dim * k + d];
1830e2ec84fSDave May }
1840e2ec84fSDave May swarm_cellid[pcnt] = e;
1850e2ec84fSDave May pcnt++;
1860e2ec84fSDave May }
1870e2ec84fSDave May }
1889566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(coor, &_coor));
189*19307e5cSMatthew G. Knepley PetscCall(DMSwarmRestoreField(dm, cellid, NULL, NULL, (void **)&swarm_cellid));
190*19307e5cSMatthew G. Knepley PetscCall(DMSwarmRestoreField(dm, coordFields[0], NULL, NULL, (void **)&swarm_coor));
1919566063dSJacob Faibussowitsch PetscCall(DMDARestoreElements(dmc, &nel, &npe, &element_list));
1920e2ec84fSDave May
1939566063dSJacob Faibussowitsch PetscCall(PetscFree(xi));
1949566063dSJacob Faibussowitsch PetscCall(PetscFree(elcoor));
19548a46eb9SPierre Jolivet for (q = 0; q < npoints_q; q++) PetscCall(PetscFree(basis[q]));
1969566063dSJacob Faibussowitsch PetscCall(PetscFree(basis));
1973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1980e2ec84fSDave May }
1990e2ec84fSDave May
private_DMSwarmInsertPointsUsingCellDM_DA(DM dm,DM celldm,DMSwarmPICLayoutType layout,PetscInt layout_param)200d71ae5a4SJacob Faibussowitsch PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_DA(DM dm, DM celldm, DMSwarmPICLayoutType layout, PetscInt layout_param)
201d71ae5a4SJacob Faibussowitsch {
2020e2ec84fSDave May DMDAElementType etype;
2030e2ec84fSDave May PetscInt dim;
2040e2ec84fSDave May
2050e2ec84fSDave May PetscFunctionBegin;
2069566063dSJacob Faibussowitsch PetscCall(DMDAGetElementType(celldm, &etype));
2079566063dSJacob Faibussowitsch PetscCall(DMGetDimension(celldm, &dim));
2080e2ec84fSDave May switch (etype) {
209d71ae5a4SJacob Faibussowitsch case DMDA_ELEMENT_P1:
210d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "DA support is not currently available for DMDA_ELEMENT_P1");
2110e2ec84fSDave May case DMDA_ELEMENT_Q1:
21208401ef6SPierre Jolivet PetscCheck(dim != 1, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Support only available for dim = 2, 3");
2139566063dSJacob Faibussowitsch PetscCall(private_DMSwarmInsertPointsUsingCellDM_DA_Q1(dm, celldm, layout_param, layout));
2140e2ec84fSDave May break;
2150e2ec84fSDave May }
2163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
2170e2ec84fSDave May }
218