Path: blob/main/src/t8_cmesh/t8_cmesh_internal/t8_cmesh_stash.c
921 views
/*1This file is part of t8code.2t8code is a C library to manage a collection (a forest) of multiple3connected adaptive space-trees of general element classes in parallel.45Copyright (C) 2015 the developers67t8code is free software; you can redistribute it and/or modify8it under the terms of the GNU General Public License as published by9the Free Software Foundation; either version 2 of the License, or10(at your option) any later version.1112t8code is distributed in the hope that it will be useful,13but WITHOUT ANY WARRANTY; without even the implied warranty of14MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the15GNU General Public License for more details.1617You should have received a copy of the GNU General Public License18along with t8code; if not, write to the Free Software Foundation, Inc.,1951 Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.20*/2122/** \file t8_cmesh_stash.c23* We define the data structures and routines for temporary storage before commit24*/2526#include <t8.h>27#include <t8_eclass.h>28#include <t8_cmesh/t8_cmesh_internal/t8_cmesh_stash.h>2930void31t8_stash_init (t8_stash_t *pstash)32{33t8_stash_t stash;34T8_ASSERT (pstash != NULL);3536stash = *pstash = T8_ALLOC (t8_stash_struct_t, 1);37sc_array_init (&stash->attributes, sizeof (t8_stash_attribute_struct_t));38sc_array_init (&stash->classes, sizeof (t8_stash_class_struct_t));39sc_array_init (&stash->joinfaces, sizeof (t8_stash_joinface_struct_t));40}4142void43t8_stash_destroy (t8_stash_t *pstash)44{45t8_stash_t stash;46t8_stash_attribute_struct_t *attr;47size_t attr_count;4849T8_ASSERT (pstash != NULL);50stash = *pstash;51sc_array_reset (&stash->classes);52sc_array_reset (&stash->joinfaces);53for (attr_count = 0; attr_count < stash->attributes.elem_count; attr_count++) {54attr = (t8_stash_attribute_struct_t *) sc_array_index (&stash->attributes, attr_count);55if (attr->is_owned) {56T8_FREE (attr->attr_data);57}58}59sc_array_reset (&stash->attributes);60T8_FREE (stash);61pstash = NULL;62}6364void65t8_stash_add_class (t8_stash_t stash, const t8_gloidx_t id, const t8_eclass_t eclass)66{67t8_stash_class_struct_t *sclass;6869T8_ASSERT (stash != NULL);70sclass = (t8_stash_class_struct_t *) sc_array_push (&stash->classes);71sclass->eclass = eclass;72sclass->id = id;73}7475/* returns -1 if treeid1 < treeid276* 0 =77* +1 >78*/79static int80t8_stash_class_compare (const void *c1, const void *c2)81{82t8_stash_class_struct_t *class1, *class2;8384class1 = (t8_stash_class_struct_t *) c1;85class2 = (t8_stash_class_struct_t *) c2;86return class1->id < class2->id ? -1 : class1->id != class2->id;87}8889void90t8_stash_class_sort (t8_stash_t stash)91{92T8_ASSERT (stash != NULL);9394sc_array_sort (&stash->classes, t8_stash_class_compare);95}9697static int98t8_stash_class_compare_index (const void *index, const void *c)99{100t8_gloidx_t index1, index2;101102index1 = *((t8_gloidx_t *) index);103index2 = ((t8_stash_class_struct_t *) c)->id;104105return index1 < index2 ? -1 : index1 != index2;106}107108ssize_t109t8_stash_class_bsearch (const t8_stash_t stash, const t8_gloidx_t tree_id)110{111return sc_array_bsearch (&stash->classes, &tree_id, t8_stash_class_compare_index);112}113114void115t8_stash_add_facejoin (t8_stash_t stash, const t8_gloidx_t gid1, const t8_gloidx_t gid2, const int face1,116const int face2, const int orientation)117{118t8_stash_joinface_struct_t *sjoin;119120T8_ASSERT (stash != NULL);121sjoin = (t8_stash_joinface_struct_t *) sc_array_push (&stash->joinfaces);122/* We insert the face connection such that join->id1 is the smaller of123* the two ids */124sjoin->face1 = gid1 <= gid2 ? face1 : face2;125sjoin->face2 = gid1 <= gid2 ? face2 : face1;126sjoin->id1 = gid1 <= gid2 ? gid1 : gid2;127sjoin->id2 = gid1 <= gid2 ? gid2 : gid1;128sjoin->orientation = orientation;129}130131/* return -1 if face1.id1 < face2.id1132* 0 =133* +1 >134*/135static int136t8_stash_facejoin_compare (const void *j1, const void *j2)137{138t8_stash_joinface_struct_t *join1, *join2;139140join1 = (t8_stash_joinface_struct_t *) j1;141join2 = (t8_stash_joinface_struct_t *) j2;142143return join1->id1 < join2->id1 ? -1 : join1->id1 != join2->id1;144}145146void147t8_stash_joinface_sort (const t8_stash_t stash)148{149T8_ASSERT (stash != NULL);150151sc_array_sort (&stash->joinfaces, t8_stash_facejoin_compare);152}153154void155t8_stash_add_attribute (const t8_stash_t stash, const t8_gloidx_t id, const int package_id, const int key,156const size_t size, void *const attr, const int copy)157{158t8_stash_attribute_struct_t *sattr;159160T8_ASSERT (stash != NULL);161sattr = (t8_stash_attribute_struct_t *) sc_array_push (&stash->attributes);162sattr->attr_size = size;163sattr->id = id;164sattr->is_owned = !copy ? 0 : 1;165sattr->key = key;166sattr->package_id = package_id;167sattr->attr_data = !copy ? attr : T8_ALLOC (char, size);168if (copy) {169memcpy (sattr->attr_data, attr, size);170}171}172173size_t174t8_stash_get_attribute_size (const t8_stash_t stash, const size_t index)175{176return ((t8_stash_attribute_struct_t *) sc_array_index (&stash->attributes, index))->attr_size;177}178179void *180t8_stash_get_attribute (const t8_stash_t stash, const size_t index)181{182return ((t8_stash_attribute_struct_t *) sc_array_index (&stash->attributes, index))->attr_data;183}184185t8_gloidx_t186t8_stash_get_attribute_tree_id (const t8_stash_t stash, const size_t index)187{188return ((t8_stash_attribute_struct_t *) sc_array_index (&stash->attributes, index))->id;189}190191int192t8_stash_get_attribute_key (const t8_stash_t stash, const size_t index)193{194return ((t8_stash_attribute_struct_t *) sc_array_index (&stash->attributes, index))->key;195}196197int198t8_stash_get_attribute_id (const t8_stash_t stash, const size_t index)199{200return ((t8_stash_attribute_struct_t *) sc_array_index (&stash->attributes, index))->package_id;201}202203int204t8_stash_attribute_is_owned (const t8_stash_t stash, const size_t index)205{206return ((t8_stash_attribute_struct_t *) sc_array_index (&stash->attributes, index))->is_owned;207}208209/* Compare two attribute entries A1 and A2.210* A1 is smaller than A2 if and only if its treeid is smaller or (if equal)211* its package id is smaller or (if also equal) its key is smaller.212*/213static int214t8_stash_attribute_compare (const void *v1, const void *v2)215{216t8_stash_attribute_struct_t *A1 = (t8_stash_attribute_struct_t *) v1;217t8_stash_attribute_struct_t *A2 = (t8_stash_attribute_struct_t *) v2;218219if (A1->id == A2->id) {220if (A1->package_id == A2->package_id) {221return A1->key < A2->key ? -1 : A1->key > A2->key;222}223return A1->package_id < A2->package_id ? -1 : 1;224}225return A1->id < A2->id ? -1 : A1->id > A2->id;226}227228/* Sort the attribute entries in the order229* (treeid, packageid, key)230*/231void232t8_stash_attribute_sort (const t8_stash_t stash)233{234sc_array_sort (&stash->attributes, t8_stash_attribute_compare);235}236237static void238t8_stash_bcast_attributes (sc_array_t *attributes, const int root, sc_MPI_Comm comm)239{240size_t num_atts, iatt, att_size, copied_bytes;241t8_stash_attribute_struct_t *att;242char *buffer;243int mpirank, mpisize, mpiret;244mpiret = sc_MPI_Comm_rank (comm, &mpirank);245SC_CHECK_MPI (mpiret);246mpiret = sc_MPI_Comm_size (comm, &mpisize);247SC_CHECK_MPI (mpiret);248249T8_ASSERT (attributes != NULL);250251num_atts = attributes->elem_count;252253/* Count the number of byte of all attributes */254att_size = 0;255for (iatt = 0; iatt < num_atts; iatt++) {256att = (t8_stash_attribute_struct_t *) sc_array_index (attributes, iatt);257att_size += att->attr_size;258}259/* Allocate buffer */260buffer = T8_ALLOC_ZERO (char, att_size);261/* Copy all attributes to send buffer */262if (mpirank == root) {263copied_bytes = 0;264for (iatt = 0; iatt < num_atts; iatt++) {265att = (t8_stash_attribute_struct_t *) sc_array_index (attributes, iatt);266memcpy (buffer + copied_bytes, att->attr_data, att->attr_size);267copied_bytes += att->attr_size;268}269}270/* broadcast buffer */271sc_MPI_Bcast (buffer, att_size, sc_MPI_BYTE, root, comm);272/* Copy attributes from buffer back to stash */273if (mpirank != root) {274copied_bytes = 0;275for (iatt = 0; iatt < num_atts; iatt++) {276att = (t8_stash_attribute_struct_t *) sc_array_index (attributes, iatt);277att->attr_data = T8_ALLOC (char, att->attr_size);278memcpy (att->attr_data, buffer + copied_bytes, att->attr_size);279copied_bytes += att->attr_size;280}281}282T8_FREE (buffer);283}284285/* bcast the data of stash on root to all procs.286* On the other procs stash_init has to be called before */287t8_stash_t288t8_stash_bcast (const t8_stash_t stash, const int root, sc_MPI_Comm comm, const size_t elem_counts[3])289{290int mpirank, mpisize, mpiret;291mpiret = sc_MPI_Comm_rank (comm, &mpirank);292SC_CHECK_MPI (mpiret);293mpiret = sc_MPI_Comm_size (comm, &mpisize);294SC_CHECK_MPI (mpiret);295296if (mpirank != root) {297sc_array_resize (&stash->attributes, elem_counts[0]);298sc_array_resize (&stash->classes, elem_counts[1]);299sc_array_resize (&stash->joinfaces, elem_counts[2]);300}301if (elem_counts[0] > 0) {302mpiret = sc_MPI_Bcast (stash->attributes.array, elem_counts[0] * sizeof (t8_stash_attribute_struct_t), sc_MPI_BYTE,3030, comm);304SC_CHECK_MPI (mpiret);305t8_stash_bcast_attributes (&stash->attributes, root, comm);306}307if (elem_counts[1] > 0) {308mpiret309= sc_MPI_Bcast (stash->classes.array, elem_counts[1] * sizeof (t8_stash_class_struct_t), sc_MPI_BYTE, 0, comm);310SC_CHECK_MPI (mpiret);311}312if (elem_counts[2] > 0) {313mpiret = sc_MPI_Bcast (stash->joinfaces.array, elem_counts[2] * sizeof (t8_stash_joinface_struct_t), sc_MPI_BYTE, 0,314comm);315SC_CHECK_MPI (mpiret);316}317return stash;318}319320int321t8_stash_is_equal (const t8_stash_t stash_a, const t8_stash_t stash_b)322{323if (stash_a == stash_b) {324return 1;325}326if (stash_a == NULL || stash_b == NULL) {327return 0;328}329return (sc_array_is_equal (&stash_a->attributes, &stash_b->attributes)330&& sc_array_is_equal (&stash_a->classes, &stash_b->classes)331&& sc_array_is_equal (&stash_a->joinfaces, &stash_b->joinfaces));332}333334335