xref: /petsc/src/dm/impls/swarm/swarmpic_da.c (revision 84467f862f2de26368b63ea79dccd665bcda30cd) !
1 #include <petscdm.h>
2 #include <petscdmda.h>                 /*I  "petscdmda.h"  I*/
3 #include <petsc/private/dmswarmimpl.h> /*I  "petscdmswarm.h"  I*/
4 #include "../src/dm/impls/swarm/data_bucket.h"
5 
6 static PetscErrorCode private_DMSwarmCreateCellLocalCoords_DA_Q1_Regular(PetscInt dim, PetscInt np[], PetscInt *_npoints, PetscReal **_xi)
7 {
8   PetscReal *xi;
9   PetscInt   d, npoints = 0, cnt;
10   PetscReal  ds[] = {0.0, 0.0, 0.0};
11   PetscInt   ii, jj, kk;
12 
13   PetscFunctionBegin;
14   switch (dim) {
15   case 1:
16     npoints = np[0];
17     break;
18   case 2:
19     npoints = np[0] * np[1];
20     break;
21   case 3:
22     npoints = np[0] * np[1] * np[2];
23     break;
24   }
25   for (d = 0; d < dim; d++) ds[d] = 2.0 / ((PetscReal)np[d]);
26 
27   PetscCall(PetscMalloc1(dim * npoints, &xi));
28   switch (dim) {
29   case 1:
30     cnt = 0;
31     for (ii = 0; ii < np[0]; ii++) {
32       xi[dim * cnt + 0] = -1.0 + 0.5 * ds[d] + ii * ds[0];
33       cnt++;
34     }
35     break;
36 
37   case 2:
38     cnt = 0;
39     for (jj = 0; jj < np[1]; jj++) {
40       for (ii = 0; ii < np[0]; ii++) {
41         xi[dim * cnt + 0] = -1.0 + 0.5 * ds[0] + ii * ds[0];
42         xi[dim * cnt + 1] = -1.0 + 0.5 * ds[1] + jj * ds[1];
43         cnt++;
44       }
45     }
46     break;
47 
48   case 3:
49     cnt = 0;
50     for (kk = 0; kk < np[2]; kk++) {
51       for (jj = 0; jj < np[1]; jj++) {
52         for (ii = 0; ii < np[0]; ii++) {
53           xi[dim * cnt + 0] = -1.0 + 0.5 * ds[0] + ii * ds[0];
54           xi[dim * cnt + 1] = -1.0 + 0.5 * ds[1] + jj * ds[1];
55           xi[dim * cnt + 2] = -1.0 + 0.5 * ds[2] + kk * ds[2];
56           cnt++;
57         }
58       }
59     }
60     break;
61   }
62   *_npoints = npoints;
63   *_xi      = xi;
64   PetscFunctionReturn(PETSC_SUCCESS);
65 }
66 
67 static PetscErrorCode private_DMSwarmCreateCellLocalCoords_DA_Q1_Gauss(PetscInt dim, PetscInt np_1d, PetscInt *_npoints, PetscReal **_xi)
68 {
69   PetscQuadrature  quadrature;
70   const PetscReal *quadrature_xi;
71   PetscReal       *xi;
72   PetscInt         d, q, npoints_q;
73 
74   PetscFunctionBegin;
75   PetscCall(PetscDTGaussTensorQuadrature(dim, 1, np_1d, -1.0, 1.0, &quadrature));
76   PetscCall(PetscQuadratureGetData(quadrature, NULL, NULL, &npoints_q, &quadrature_xi, NULL));
77   PetscCall(PetscMalloc1(dim * npoints_q, &xi));
78   for (q = 0; q < npoints_q; q++) {
79     for (d = 0; d < dim; d++) xi[dim * q + d] = quadrature_xi[dim * q + d];
80   }
81   PetscCall(PetscQuadratureDestroy(&quadrature));
82   *_npoints = npoints_q;
83   *_xi      = xi;
84   PetscFunctionReturn(PETSC_SUCCESS);
85 }
86 
87 static PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_DA_Q1(DM dm, DM dmc, PetscInt npoints, DMSwarmPICLayoutType layout)
88 {
89   PetscInt           dim, npoints_q;
90   PetscInt           nel, npe, e, q, k, d;
91   const PetscInt    *element_list;
92   PetscReal        **basis;
93   PetscReal         *xi;
94   Vec                coor;
95   const PetscScalar *_coor;
96   PetscReal         *elcoor;
97   PetscReal         *swarm_coor;
98   PetscInt          *swarm_cellid;
99   PetscInt           pcnt;
100   const char        *coordname;
101 
102   PetscFunctionBegin;
103   PetscCall(DMGetDimension(dm, &dim));
104   switch (layout) {
105   case DMSWARMPIC_LAYOUT_REGULAR: {
106     PetscInt np_dir[3];
107     np_dir[0] = np_dir[1] = np_dir[2] = npoints;
108     PetscCall(private_DMSwarmCreateCellLocalCoords_DA_Q1_Regular(dim, np_dir, &npoints_q, &xi));
109   } break;
110   case DMSWARMPIC_LAYOUT_GAUSS:
111     PetscCall(private_DMSwarmCreateCellLocalCoords_DA_Q1_Gauss(dim, npoints, &npoints_q, &xi));
112     break;
113 
114   case DMSWARMPIC_LAYOUT_SUBDIVISION: {
115     PetscInt s, nsub;
116     PetscInt np_dir[3];
117     nsub      = npoints;
118     np_dir[0] = 1;
119     for (s = 0; s < nsub; s++) np_dir[0] *= 2;
120     np_dir[1] = np_dir[0];
121     np_dir[2] = np_dir[0];
122     PetscCall(private_DMSwarmCreateCellLocalCoords_DA_Q1_Regular(dim, np_dir, &npoints_q, &xi));
123   } break;
124   default:
125     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "A valid DMSwarmPIC layout must be provided");
126   }
127 
128   PetscCall(DMDAGetElements(dmc, &nel, &npe, &element_list));
129   PetscCall(PetscMalloc1(dim * npe, &elcoor));
130   PetscCall(PetscMalloc1(npoints_q, &basis));
131   for (q = 0; q < npoints_q; q++) {
132     PetscCall(PetscMalloc1(npe, &basis[q]));
133 
134     switch (dim) {
135     case 1:
136       basis[q][0] = 0.5 * (1.0 - xi[dim * q + 0]);
137       basis[q][1] = 0.5 * (1.0 + xi[dim * q + 0]);
138       break;
139     case 2:
140       basis[q][0] = 0.25 * (1.0 - xi[dim * q + 0]) * (1.0 - xi[dim * q + 1]);
141       basis[q][1] = 0.25 * (1.0 + xi[dim * q + 0]) * (1.0 - xi[dim * q + 1]);
142       basis[q][2] = 0.25 * (1.0 + xi[dim * q + 0]) * (1.0 + xi[dim * q + 1]);
143       basis[q][3] = 0.25 * (1.0 - xi[dim * q + 0]) * (1.0 + xi[dim * q + 1]);
144       break;
145 
146     case 3:
147       basis[q][0] = 0.125 * (1.0 - xi[dim * q + 0]) * (1.0 - xi[dim * q + 1]) * (1.0 - xi[dim * q + 2]);
148       basis[q][1] = 0.125 * (1.0 + xi[dim * q + 0]) * (1.0 - xi[dim * q + 1]) * (1.0 - xi[dim * q + 2]);
149       basis[q][2] = 0.125 * (1.0 + xi[dim * q + 0]) * (1.0 + xi[dim * q + 1]) * (1.0 - xi[dim * q + 2]);
150       basis[q][3] = 0.125 * (1.0 - xi[dim * q + 0]) * (1.0 + xi[dim * q + 1]) * (1.0 - xi[dim * q + 2]);
151       basis[q][4] = 0.125 * (1.0 - xi[dim * q + 0]) * (1.0 - xi[dim * q + 1]) * (1.0 + xi[dim * q + 2]);
152       basis[q][5] = 0.125 * (1.0 + xi[dim * q + 0]) * (1.0 - xi[dim * q + 1]) * (1.0 + xi[dim * q + 2]);
153       basis[q][6] = 0.125 * (1.0 + xi[dim * q + 0]) * (1.0 + xi[dim * q + 1]) * (1.0 + xi[dim * q + 2]);
154       basis[q][7] = 0.125 * (1.0 - xi[dim * q + 0]) * (1.0 + xi[dim * q + 1]) * (1.0 + xi[dim * q + 2]);
155       break;
156     }
157   }
158 
159   PetscCall(DMSwarmSetLocalSizes(dm, npoints_q * nel, -1));
160   PetscCall(DMSwarmGetCoordinateField(dm, &coordname));
161   PetscCall(DMSwarmGetField(dm, coordname, NULL, NULL, (void **)&swarm_coor));
162   PetscCall(DMSwarmGetField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid));
163 
164   PetscCall(DMGetCoordinatesLocal(dmc, &coor));
165   PetscCall(VecGetArrayRead(coor, &_coor));
166   pcnt = 0;
167   for (e = 0; e < nel; e++) {
168     const PetscInt *element = &element_list[npe * e];
169 
170     for (k = 0; k < npe; k++) {
171       for (d = 0; d < dim; d++) elcoor[dim * k + d] = PetscRealPart(_coor[dim * element[k] + d]);
172     }
173 
174     for (q = 0; q < npoints_q; q++) {
175       for (d = 0; d < dim; d++) swarm_coor[dim * pcnt + d] = 0.0;
176       for (k = 0; k < npe; k++) {
177         for (d = 0; d < dim; d++) swarm_coor[dim * pcnt + d] += basis[q][k] * elcoor[dim * k + d];
178       }
179       swarm_cellid[pcnt] = e;
180       pcnt++;
181     }
182   }
183   PetscCall(VecRestoreArrayRead(coor, &_coor));
184   PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid));
185   PetscCall(DMSwarmRestoreField(dm, coordname, NULL, NULL, (void **)&swarm_coor));
186   PetscCall(DMDARestoreElements(dmc, &nel, &npe, &element_list));
187 
188   PetscCall(PetscFree(xi));
189   PetscCall(PetscFree(elcoor));
190   for (q = 0; q < npoints_q; q++) PetscCall(PetscFree(basis[q]));
191   PetscCall(PetscFree(basis));
192   PetscFunctionReturn(PETSC_SUCCESS);
193 }
194 
195 PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_DA(DM dm, DM celldm, DMSwarmPICLayoutType layout, PetscInt layout_param)
196 {
197   DMDAElementType etype;
198   PetscInt        dim;
199 
200   PetscFunctionBegin;
201   PetscCall(DMDAGetElementType(celldm, &etype));
202   PetscCall(DMGetDimension(celldm, &dim));
203   switch (etype) {
204   case DMDA_ELEMENT_P1:
205     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "DA support is not currently available for DMDA_ELEMENT_P1");
206   case DMDA_ELEMENT_Q1:
207     PetscCheck(dim != 1, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Support only available for dim = 2, 3");
208     PetscCall(private_DMSwarmInsertPointsUsingCellDM_DA_Q1(dm, celldm, layout_param, layout));
209     break;
210   }
211   PetscFunctionReturn(PETSC_SUCCESS);
212 }
213