Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
DLR-AMR
GitHub Repository: DLR-AMR/t8code
Path: blob/main/src/t8_cmesh/t8_cmesh_internal/t8_cmesh_stash.c
921 views
1
/*
2
This file is part of t8code.
3
t8code is a C library to manage a collection (a forest) of multiple
4
connected adaptive space-trees of general element classes in parallel.
5
6
Copyright (C) 2015 the developers
7
8
t8code is free software; you can redistribute it and/or modify
9
it under the terms of the GNU General Public License as published by
10
the Free Software Foundation; either version 2 of the License, or
11
(at your option) any later version.
12
13
t8code is distributed in the hope that it will be useful,
14
but WITHOUT ANY WARRANTY; without even the implied warranty of
15
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16
GNU General Public License for more details.
17
18
You should have received a copy of the GNU General Public License
19
along with t8code; if not, write to the Free Software Foundation, Inc.,
20
51 Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
21
*/
22
23
/** \file t8_cmesh_stash.c
24
* We define the data structures and routines for temporary storage before commit
25
*/
26
27
#include <t8.h>
28
#include <t8_eclass.h>
29
#include <t8_cmesh/t8_cmesh_internal/t8_cmesh_stash.h>
30
31
void
32
t8_stash_init (t8_stash_t *pstash)
33
{
34
t8_stash_t stash;
35
T8_ASSERT (pstash != NULL);
36
37
stash = *pstash = T8_ALLOC (t8_stash_struct_t, 1);
38
sc_array_init (&stash->attributes, sizeof (t8_stash_attribute_struct_t));
39
sc_array_init (&stash->classes, sizeof (t8_stash_class_struct_t));
40
sc_array_init (&stash->joinfaces, sizeof (t8_stash_joinface_struct_t));
41
}
42
43
void
44
t8_stash_destroy (t8_stash_t *pstash)
45
{
46
t8_stash_t stash;
47
t8_stash_attribute_struct_t *attr;
48
size_t attr_count;
49
50
T8_ASSERT (pstash != NULL);
51
stash = *pstash;
52
sc_array_reset (&stash->classes);
53
sc_array_reset (&stash->joinfaces);
54
for (attr_count = 0; attr_count < stash->attributes.elem_count; attr_count++) {
55
attr = (t8_stash_attribute_struct_t *) sc_array_index (&stash->attributes, attr_count);
56
if (attr->is_owned) {
57
T8_FREE (attr->attr_data);
58
}
59
}
60
sc_array_reset (&stash->attributes);
61
T8_FREE (stash);
62
pstash = NULL;
63
}
64
65
void
66
t8_stash_add_class (t8_stash_t stash, const t8_gloidx_t id, const t8_eclass_t eclass)
67
{
68
t8_stash_class_struct_t *sclass;
69
70
T8_ASSERT (stash != NULL);
71
sclass = (t8_stash_class_struct_t *) sc_array_push (&stash->classes);
72
sclass->eclass = eclass;
73
sclass->id = id;
74
}
75
76
/* returns -1 if treeid1 < treeid2
77
* 0 =
78
* +1 >
79
*/
80
static int
81
t8_stash_class_compare (const void *c1, const void *c2)
82
{
83
t8_stash_class_struct_t *class1, *class2;
84
85
class1 = (t8_stash_class_struct_t *) c1;
86
class2 = (t8_stash_class_struct_t *) c2;
87
return class1->id < class2->id ? -1 : class1->id != class2->id;
88
}
89
90
void
91
t8_stash_class_sort (t8_stash_t stash)
92
{
93
T8_ASSERT (stash != NULL);
94
95
sc_array_sort (&stash->classes, t8_stash_class_compare);
96
}
97
98
static int
99
t8_stash_class_compare_index (const void *index, const void *c)
100
{
101
t8_gloidx_t index1, index2;
102
103
index1 = *((t8_gloidx_t *) index);
104
index2 = ((t8_stash_class_struct_t *) c)->id;
105
106
return index1 < index2 ? -1 : index1 != index2;
107
}
108
109
ssize_t
110
t8_stash_class_bsearch (const t8_stash_t stash, const t8_gloidx_t tree_id)
111
{
112
return sc_array_bsearch (&stash->classes, &tree_id, t8_stash_class_compare_index);
113
}
114
115
void
116
t8_stash_add_facejoin (t8_stash_t stash, const t8_gloidx_t gid1, const t8_gloidx_t gid2, const int face1,
117
const int face2, const int orientation)
118
{
119
t8_stash_joinface_struct_t *sjoin;
120
121
T8_ASSERT (stash != NULL);
122
sjoin = (t8_stash_joinface_struct_t *) sc_array_push (&stash->joinfaces);
123
/* We insert the face connection such that join->id1 is the smaller of
124
* the two ids */
125
sjoin->face1 = gid1 <= gid2 ? face1 : face2;
126
sjoin->face2 = gid1 <= gid2 ? face2 : face1;
127
sjoin->id1 = gid1 <= gid2 ? gid1 : gid2;
128
sjoin->id2 = gid1 <= gid2 ? gid2 : gid1;
129
sjoin->orientation = orientation;
130
}
131
132
/* return -1 if face1.id1 < face2.id1
133
* 0 =
134
* +1 >
135
*/
136
static int
137
t8_stash_facejoin_compare (const void *j1, const void *j2)
138
{
139
t8_stash_joinface_struct_t *join1, *join2;
140
141
join1 = (t8_stash_joinface_struct_t *) j1;
142
join2 = (t8_stash_joinface_struct_t *) j2;
143
144
return join1->id1 < join2->id1 ? -1 : join1->id1 != join2->id1;
145
}
146
147
void
148
t8_stash_joinface_sort (const t8_stash_t stash)
149
{
150
T8_ASSERT (stash != NULL);
151
152
sc_array_sort (&stash->joinfaces, t8_stash_facejoin_compare);
153
}
154
155
void
156
t8_stash_add_attribute (const t8_stash_t stash, const t8_gloidx_t id, const int package_id, const int key,
157
const size_t size, void *const attr, const int copy)
158
{
159
t8_stash_attribute_struct_t *sattr;
160
161
T8_ASSERT (stash != NULL);
162
sattr = (t8_stash_attribute_struct_t *) sc_array_push (&stash->attributes);
163
sattr->attr_size = size;
164
sattr->id = id;
165
sattr->is_owned = !copy ? 0 : 1;
166
sattr->key = key;
167
sattr->package_id = package_id;
168
sattr->attr_data = !copy ? attr : T8_ALLOC (char, size);
169
if (copy) {
170
memcpy (sattr->attr_data, attr, size);
171
}
172
}
173
174
size_t
175
t8_stash_get_attribute_size (const t8_stash_t stash, const size_t index)
176
{
177
return ((t8_stash_attribute_struct_t *) sc_array_index (&stash->attributes, index))->attr_size;
178
}
179
180
void *
181
t8_stash_get_attribute (const t8_stash_t stash, const size_t index)
182
{
183
return ((t8_stash_attribute_struct_t *) sc_array_index (&stash->attributes, index))->attr_data;
184
}
185
186
t8_gloidx_t
187
t8_stash_get_attribute_tree_id (const t8_stash_t stash, const size_t index)
188
{
189
return ((t8_stash_attribute_struct_t *) sc_array_index (&stash->attributes, index))->id;
190
}
191
192
int
193
t8_stash_get_attribute_key (const t8_stash_t stash, const size_t index)
194
{
195
return ((t8_stash_attribute_struct_t *) sc_array_index (&stash->attributes, index))->key;
196
}
197
198
int
199
t8_stash_get_attribute_id (const t8_stash_t stash, const size_t index)
200
{
201
return ((t8_stash_attribute_struct_t *) sc_array_index (&stash->attributes, index))->package_id;
202
}
203
204
int
205
t8_stash_attribute_is_owned (const t8_stash_t stash, const size_t index)
206
{
207
return ((t8_stash_attribute_struct_t *) sc_array_index (&stash->attributes, index))->is_owned;
208
}
209
210
/* Compare two attribute entries A1 and A2.
211
* A1 is smaller than A2 if and only if its treeid is smaller or (if equal)
212
* its package id is smaller or (if also equal) its key is smaller.
213
*/
214
static int
215
t8_stash_attribute_compare (const void *v1, const void *v2)
216
{
217
t8_stash_attribute_struct_t *A1 = (t8_stash_attribute_struct_t *) v1;
218
t8_stash_attribute_struct_t *A2 = (t8_stash_attribute_struct_t *) v2;
219
220
if (A1->id == A2->id) {
221
if (A1->package_id == A2->package_id) {
222
return A1->key < A2->key ? -1 : A1->key > A2->key;
223
}
224
return A1->package_id < A2->package_id ? -1 : 1;
225
}
226
return A1->id < A2->id ? -1 : A1->id > A2->id;
227
}
228
229
/* Sort the attribute entries in the order
230
* (treeid, packageid, key)
231
*/
232
void
233
t8_stash_attribute_sort (const t8_stash_t stash)
234
{
235
sc_array_sort (&stash->attributes, t8_stash_attribute_compare);
236
}
237
238
static void
239
t8_stash_bcast_attributes (sc_array_t *attributes, const int root, sc_MPI_Comm comm)
240
{
241
size_t num_atts, iatt, att_size, copied_bytes;
242
t8_stash_attribute_struct_t *att;
243
char *buffer;
244
int mpirank, mpisize, mpiret;
245
mpiret = sc_MPI_Comm_rank (comm, &mpirank);
246
SC_CHECK_MPI (mpiret);
247
mpiret = sc_MPI_Comm_size (comm, &mpisize);
248
SC_CHECK_MPI (mpiret);
249
250
T8_ASSERT (attributes != NULL);
251
252
num_atts = attributes->elem_count;
253
254
/* Count the number of byte of all attributes */
255
att_size = 0;
256
for (iatt = 0; iatt < num_atts; iatt++) {
257
att = (t8_stash_attribute_struct_t *) sc_array_index (attributes, iatt);
258
att_size += att->attr_size;
259
}
260
/* Allocate buffer */
261
buffer = T8_ALLOC_ZERO (char, att_size);
262
/* Copy all attributes to send buffer */
263
if (mpirank == root) {
264
copied_bytes = 0;
265
for (iatt = 0; iatt < num_atts; iatt++) {
266
att = (t8_stash_attribute_struct_t *) sc_array_index (attributes, iatt);
267
memcpy (buffer + copied_bytes, att->attr_data, att->attr_size);
268
copied_bytes += att->attr_size;
269
}
270
}
271
/* broadcast buffer */
272
sc_MPI_Bcast (buffer, att_size, sc_MPI_BYTE, root, comm);
273
/* Copy attributes from buffer back to stash */
274
if (mpirank != root) {
275
copied_bytes = 0;
276
for (iatt = 0; iatt < num_atts; iatt++) {
277
att = (t8_stash_attribute_struct_t *) sc_array_index (attributes, iatt);
278
att->attr_data = T8_ALLOC (char, att->attr_size);
279
memcpy (att->attr_data, buffer + copied_bytes, att->attr_size);
280
copied_bytes += att->attr_size;
281
}
282
}
283
T8_FREE (buffer);
284
}
285
286
/* bcast the data of stash on root to all procs.
287
* On the other procs stash_init has to be called before */
288
t8_stash_t
289
t8_stash_bcast (const t8_stash_t stash, const int root, sc_MPI_Comm comm, const size_t elem_counts[3])
290
{
291
int mpirank, mpisize, mpiret;
292
mpiret = sc_MPI_Comm_rank (comm, &mpirank);
293
SC_CHECK_MPI (mpiret);
294
mpiret = sc_MPI_Comm_size (comm, &mpisize);
295
SC_CHECK_MPI (mpiret);
296
297
if (mpirank != root) {
298
sc_array_resize (&stash->attributes, elem_counts[0]);
299
sc_array_resize (&stash->classes, elem_counts[1]);
300
sc_array_resize (&stash->joinfaces, elem_counts[2]);
301
}
302
if (elem_counts[0] > 0) {
303
mpiret = sc_MPI_Bcast (stash->attributes.array, elem_counts[0] * sizeof (t8_stash_attribute_struct_t), sc_MPI_BYTE,
304
0, comm);
305
SC_CHECK_MPI (mpiret);
306
t8_stash_bcast_attributes (&stash->attributes, root, comm);
307
}
308
if (elem_counts[1] > 0) {
309
mpiret
310
= sc_MPI_Bcast (stash->classes.array, elem_counts[1] * sizeof (t8_stash_class_struct_t), sc_MPI_BYTE, 0, comm);
311
SC_CHECK_MPI (mpiret);
312
}
313
if (elem_counts[2] > 0) {
314
mpiret = sc_MPI_Bcast (stash->joinfaces.array, elem_counts[2] * sizeof (t8_stash_joinface_struct_t), sc_MPI_BYTE, 0,
315
comm);
316
SC_CHECK_MPI (mpiret);
317
}
318
return stash;
319
}
320
321
int
322
t8_stash_is_equal (const t8_stash_t stash_a, const t8_stash_t stash_b)
323
{
324
if (stash_a == stash_b) {
325
return 1;
326
}
327
if (stash_a == NULL || stash_b == NULL) {
328
return 0;
329
}
330
return (sc_array_is_equal (&stash_a->attributes, &stash_b->attributes)
331
&& sc_array_is_equal (&stash_a->classes, &stash_b->classes)
332
&& sc_array_is_equal (&stash_a->joinfaces, &stash_b->joinfaces));
333
}
334
335