62 #pragma GCC diagnostic ignored "-Wold-style-cast" 90 ZOLTAN_ID_PTR globalID, ZOLTAN_ID_PTR localID,
91 int wgt_dim,
float *obj_wgts,
int *ierr)
104 for (Foam::label i=0; i<
mesh.nCells(); i++)
114 ZOLTAN_ID_PTR globalID, ZOLTAN_ID_PTR localID,
115 int *numEdges,
int *ierr)
119 if ((sizeGID != 1) || (sizeLID != 1) || (num_obj !=
mesh.nCells()))
121 *ierr = ZOLTAN_FATAL;
125 for (Foam::label i=0; i < num_obj ;i++)
127 Foam::label celli = localID[i];
132 if (
mesh.isInternalFace(cFaces[cFacei]))
144 static void get_edge_list(
void *data,
int sizeGID,
int sizeLID,
145 int num_obj, ZOLTAN_ID_PTR globalID, ZOLTAN_ID_PTR localID,
147 ZOLTAN_ID_PTR nborGID,
int *nborProc,
148 int wgt_dim,
float *ewgts,
int *ierr)
156 || (num_obj !=
mesh.nCells())
160 *ierr = ZOLTAN_FATAL;
164 ZOLTAN_ID_TYPE* nextNbor = nborGID;
165 int* nextProc = nborProc;
166 float* nextWgt = ewgts;
168 for (Foam::label i=0; i < num_obj; i++)
170 Foam::label celli = localID[i];
177 Foam::label facei = cFaces[cFacei];
178 if (
mesh.isInternalFace(facei))
180 Foam::label nbr =
mesh.faceOwner()[facei];
183 nbr =
mesh.faceNeighbour()[facei];
193 if (
n != num_edges[i])
195 *ierr = ZOLTAN_FATAL;
218 ZOLTAN_ID_PTR global_ids,
219 ZOLTAN_ID_PTR local_ids,
229 (num_gid_entries != 1)
230 || (num_lid_entries != 1)
231 || (num_obj !=
mesh.nCells())
232 || (num_dim !=
mesh.nSolutionD())
235 *ierr = ZOLTAN_FATAL;
239 double*
p = geom_vec;
246 for (Foam::label celli = 0; celli < num_obj; celli++)
250 for (
Foam::direction cmpt = 0; cmpt < Foam::vector::nComponents; cmpt++)
267 coeffsDict_(
dict.optionalSubDict(typeName+
"Coeffs"))
280 args[0] =
"zoltanRenumber";
284 for (label i = 0; i < argc; i++)
286 argv[i] = strdup(
args[i].c_str());
290 int rc = Zoltan_Initialize(argc, argv, &ver);
300 struct Zoltan_Struct *zz = Zoltan_Create(MPI_COMM_WORLD);
305 for (
const entry& dEntry : coeffsDict_)
307 if (!dEntry.isDict())
309 const word&
key = dEntry.keyword();
310 const word value(dEntry.get<
word>());
312 Info<< typeName <<
" : setting parameter " <<
key 313 <<
" to " << value <<
endl;
315 Zoltan_Set_Param(zz,
key.c_str(), value.c_str());
340 wantedCells[i] = globalCells.toGlobal(i);
345 int err = Zoltan_Order
354 if (err != ZOLTAN_OK)
361 for (label i = 0; i < argc; i++)
370 order[i] = oldToNew[i];
static int get_mesh_dim(void *data, int *ierr)
errorManipArg< error, int > exit(error &err, const int errNo=1)
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Ostream & endl(Ostream &os)
Add newline and flush stream.
static void get_geom_list(void *data, int num_gid_entries, int num_lid_entries, int num_obj, ZOLTAN_ID_PTR global_ids, ZOLTAN_ID_PTR local_ids, int num_dim, double *geom_vec, int *ierr)
virtual labelList renumber(const pointField &) const
Return the order in which cells need to be visited (ie. from ordered back to original cell label)...
Macros for easy insertion into run-time selection tables.
#define forAll(list, i)
Loop across all elements in list.
Abstract base class for renumbering.
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
static void get_edge_list(void *data, int sizeGID, int sizeLID, int num_obj, ZOLTAN_ID_PTR globalID, ZOLTAN_ID_PTR localID, int *num_edges, ZOLTAN_ID_PTR nborGID, int *nborProc, int wgt_dim, float *ewgts, int *ierr)
A class for handling words, derived from Foam::string.
label size() const noexcept
The number of arguments.
const globalMeshData & globalData() const
Return parallel info.
Templated 3D Vector derived from VectorSpace adding construction from 3 components, element access using x(), y() and z() member functions and the inner-product (dot-product) and cross-product operators.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
defineTypeNameAndDebug(combustionModel, 0)
label nTotalCells() const noexcept
Return total number of cells in decomposed mesh.
static int get_number_of_vertices(void *data, int *ierr)
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
static void get_vertex_list(void *data, int sizeGID, int sizeLID, ZOLTAN_ID_PTR globalID, ZOLTAN_ID_PTR localID, int wgt_dim, float *obj_wgts, int *ierr)
label nCells() const noexcept
Number of mesh cells.
auto key(const Type &t) -> typename std::enable_if< std::is_enum< Type >::value, typename std::underlying_type< Type >::type >::type
A cell is defined as a list of faces with extra functionality.
messageStream Info
Information stream (stdout output on master, null elsewhere)
Mesh consisting of general polyhedral cells.
Foam::argList args(argc, argv)
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
A keyword and a list of tokens is an 'entry'.
static void get_num_edges_list(void *data, int sizeGID, int sizeLID, int num_obj, ZOLTAN_ID_PTR globalID, ZOLTAN_ID_PTR localID, int *numEdges, int *ierr)