Skip to content

Commit

Permalink
Merge pull request #89 from imadhammani/main
Browse files Browse the repository at this point in the history
XCore: AdaptMesh OK without MPI
  • Loading branch information
benoit128 authored Jul 18, 2024
2 parents 66d4acf + 36d48f0 commit 95eb4ec
Show file tree
Hide file tree
Showing 18 changed files with 93 additions and 80 deletions.
8 changes: 3 additions & 5 deletions Cassiopee/XCore/XCore/AdaptMesh/AdaptMesh_Adapt.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -45,12 +45,14 @@ void cache_prerefinement_data(Mesh *M)
static
void print_postrefinement_data(Mesh *M)
{
Int gnc = 0;
Int gnc = M->nc;
MPI_Allreduce(&M->nc, &gnc, 1, XMPI_INT, MPI_SUM, MPI_COMM_WORLD);
if (M->pid == 0) {
printf(" Total cells after refinement: " SF_D_ "\n", gnc);
}

if (M->npc == 1) return;

Float balanced = gnc / (Float) M->npc;
Float my_imbalance = fabs((M->nc - balanced) / (Float)balanced * 100.0);
Float max_imbalance;
Expand Down Expand Up @@ -120,12 +122,8 @@ PyObject *K_XCORE::AdaptMesh_Adapt(PyObject *self, PyObject *args)

if (M->pid == 0) puts(" Conformizing face edges...");
Mesh_conformize_face_edge(M);

//if (M->pid == 0) puts(" Exporting CGNS array...");
//PyObject *karray = Mesh_export_karray(M);

if (M->pid == 0) puts(" Done.");

return Py_None;
//return karray;
}
8 changes: 6 additions & 2 deletions Cassiopee/XCore/XCore/AdaptMesh/AdaptMesh_Init.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -103,13 +103,17 @@ PyObject *K_XCORE::AdaptMesh_Init(PyObject *self, PyObject *args)

M->bps[i].gid = PyLong_AsLong(PyList_GetItem(PATCH, 1));

#if PY_VERSION_HEX >= 0x03000000
const char *bcname = PyUnicode_AsUTF8(PyList_GetItem(PATCH, 2));
const char *bctype = PyUnicode_AsUTF8(PyList_GetItem(PATCH, 3));
#else
const char *bcname = PyString_AsString(PyList_GetItem(PATCH, 2));
const char *bctype = PyString_AsString(PyList_GetItem(PATCH, 3));
#endif

M->bps[i].name = (char *)XMALLOC(strlen(bcname)+1);
strcpy(M->bps[i].name, bcname);

const char *bctype = PyUnicode_AsUTF8(PyList_GetItem(PATCH, 3));

M->bps[i].type = (char *)XMALLOC(strlen(bctype)+1);
strcpy(M->bps[i].type, bctype);

Expand Down
4 changes: 3 additions & 1 deletion Cassiopee/XCore/XCore/AdaptMesh/AdaptMesh_LoadBalance.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -34,6 +34,8 @@ PyObject *K_XCORE::AdaptMesh_LoadBalance(PyObject *self, PyObject *args)

Mesh *M = (Mesh *)PyCapsule_GetPointer(MESH, "AdaptMesh");

if (M->npc == 1) return Py_None;

// Make global cell-cell connectivity
Mesh_make_cell_cells(M);

Expand All @@ -45,7 +47,7 @@ PyObject *K_XCORE::AdaptMesh_LoadBalance(PyObject *self, PyObject *args)
return NULL;
}

if (M->pid == 0) puts("OK LoadBalance");
if (M->pid == 0) puts("OK LoadBalance.");

return Py_None;
}
2 changes: 1 addition & 1 deletion Cassiopee/XCore/XCore/AdaptMesh/Karray.h
Original file line number Diff line number Diff line change
Expand Up @@ -19,7 +19,7 @@
#pragma once

#include "xcore.h"
#include "../common/common.h"
#include "common/common.h"

struct Karray {
// Reference to the python array
Expand Down
18 changes: 11 additions & 7 deletions Cassiopee/XCore/XCore/AdaptMesh/Mesh.h
Original file line number Diff line number Diff line change
Expand Up @@ -18,12 +18,12 @@
*/
#pragma once

#include <mpi.h>
#include "Mpi.h"
#include <map>
#include <array>

#include "xcore.h"
#include "../common/common.h"
#include "common/common.h"

#define HEXA 0
#define TETRA 1
Expand Down Expand Up @@ -152,9 +152,9 @@ struct Mesh {

/* Parallel */

Int pid;
Int npc;
Int nrq;
int pid;
int npc;
int nrq;
MPI_Request *reqs;

Int *l2gc;
Expand Down Expand Up @@ -331,11 +331,15 @@ PyObject *Mesh_export_karray(Mesh *M, int conformize);

/* Parallel */

void Mesh_comm_waitall(Mesh *M);
inline
void Mesh_comm_waitall(Mesh *M)
{
MPI_Waitall(M->nrq, M->reqs, MPI_STATUSES_IGNORE);
M->nrq = 0;
}

Int Mesh_load_balance(Mesh *M);


/* Adaptation */

Int Mesh_smooth_cref(Mesh *M);
Expand Down
24 changes: 2 additions & 22 deletions Cassiopee/XCore/XCore/AdaptMesh/MeshComm.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -20,8 +20,8 @@
#include <limits>

#include "Mesh.h"
#include "common/mem.h"
#include "scotch/ptscotch.h"
#include "../common/mem.h"

#define TOL 1e-12

Expand Down Expand Up @@ -94,20 +94,6 @@ Int *compute_cell_weights(Mesh *M)
cwgts[i] = ((weights[i]/min_weight - 1)*range_scale) + 1;
}

/*
for (Int i = 0; i < M->nc; i++) {
Int cval = M->cref[i];
cwgts[i] = 1;
if (cval > 0) cwgts[i] += gnc;
//if (cval == 0) cwgts[i] = 1;
//else cwgts[i] = 8 + 8*20;
//cwgts[i] = 1 + ((1 << (3*cval)) - 1);
}
*/

return cwgts;
}

Expand Down Expand Up @@ -1546,10 +1532,4 @@ Int Mesh_load_balance(Mesh *M)
MPI_Allreduce(&M->nf, &gnf, 1, XMPI_INT, MPI_SUM, MPI_COMM_WORLD);

return gret;
}

void Mesh_comm_waitall(Mesh *M)
{
MPI_Waitall(M->nrq, M->reqs, MPI_STATUSES_IGNORE);
M->nrq = 0;
}
}
2 changes: 2 additions & 0 deletions Cassiopee/XCore/XCore/AdaptMesh/MeshConformize.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -70,6 +70,8 @@ void Mesh_conformize_face_edge(Mesh *M)
}
}

if (M->npc == 1) return;

// Exchange franges at the interface

for (Int i = 0; i < M->npp; i++) {
Expand Down
13 changes: 7 additions & 6 deletions Cassiopee/XCore/XCore/AdaptMesh/MeshConnectivity.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -168,12 +168,13 @@ void Mesh_make_edge_connectivity(Mesh *M)
}
}


void Mesh_update_global_cell_ids(Mesh *M)
{
if (M->npc == 1) return;

Int new_cells = M->nc - M->nc_old;

Int first_new_cell = 0;
Int first_new_cell = new_cells;
MPI_Scan(&new_cells, &first_new_cell, 1, XMPI_INT, MPI_SUM, MPI_COMM_WORLD);

Int gnc_old = 0;
Expand All @@ -196,8 +197,6 @@ void Mesh_update_global_cell_ids(Mesh *M)
M->g2lc[gcid] = lcid;
}
}

//MPI_Barrier(MPI_COMM_WORLD);
}

void Mesh_update_ppatches(Mesh *M)
Expand Down Expand Up @@ -365,6 +364,8 @@ void Mesh_update_bpatches(Mesh *M)

void Mesh_update_global_face_ids(Mesh *M)
{
if (M->npc == 1) return;

M->l2gf = (Int *)XRESIZE(M->l2gf, M->nf * sizeof(Int));
for (Int i = M->nf_old; i < M->nf; i++) M->l2gf[i] = -1;

Expand All @@ -377,10 +378,10 @@ void Mesh_update_global_face_ids(Mesh *M)

Int ndup = M->nf - nfree;

Int gnfree;
Int gnfree = nfree;
MPI_Allreduce(&nfree, &gnfree, 1, XMPI_INT, MPI_SUM, MPI_COMM_WORLD);

Int gndup;
Int gndup = ndup;
MPI_Allreduce(&ndup, &gndup, 1, XMPI_INT, MPI_SUM, MPI_COMM_WORLD);

gndup /= 2;
Expand Down
8 changes: 4 additions & 4 deletions Cassiopee/XCore/XCore/AdaptMesh/MeshInit.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -68,10 +68,10 @@ Mesh::Mesh()
nf_old = 0;

/* Parallel */
int pid_; int npc_;
MPI_Comm_rank(MPI_COMM_WORLD, &pid_);
MPI_Comm_size(MPI_COMM_WORLD, &npc_);
pid = E_Int(pid_); npc = E_Int(npc_);
pid = 0;
npc = 1;
MPI_Comm_rank(MPI_COMM_WORLD, &pid);
MPI_Comm_size(MPI_COMM_WORLD, &npc);

nrq = 0;
reqs = (MPI_Request *)XMALLOC(2 * npc * sizeof(MPI_Request));
Expand Down
2 changes: 1 addition & 1 deletion Cassiopee/XCore/XCore/AdaptMesh/MeshOrient.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -21,7 +21,7 @@
#include <stack>

#include "Mesh.h"
#include "../common/mem.h"
#include "common/mem.h"
#include "Edge.h"

#define INTERNAL 0
Expand Down
4 changes: 2 additions & 2 deletions Cassiopee/XCore/XCore/AdaptMesh/MeshSmooth.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -19,7 +19,7 @@
#include <stack>

#include "Mesh.h"
#include "../common/mem.h"
#include "common/mem.h"

inline
Int Mesh_get_cnei(Mesh *M, Int cid, Int fid)
Expand Down Expand Up @@ -125,7 +125,7 @@ Int Mesh_smooth_cref(Mesh *M)
}
}

Int gstop;
Int gstop = 0;
MPI_Allreduce(&lstop, &gstop, 1, XMPI_INT, MPI_SUM, MPI_COMM_WORLD);

if (gstop == 0) break;
Expand Down
24 changes: 24 additions & 0 deletions Cassiopee/XCore/XCore/AdaptMesh/Mpi.h
Original file line number Diff line number Diff line change
@@ -0,0 +1,24 @@
#pragma once

#ifdef _MPI

#include <mpi.h>

#else

#define MPI_Comm_rank(a, b)
#define MPI_Comm_size(a, b)

#define MPI_Request char

#define MPI_Isend(a, b, c, d, e, f, g)
#define MPI_Irecv(a, b, c, d, e, f, g)
#define MPI_Waitall(a)

#define MPI_Barrier(a)
#define MPI_Allreduce(a, b, c, d, e, f)
#define MPI_Waitall(a, b, c)

#define MPI_Scan(a, b, c, d, e, f)

#endif
2 changes: 1 addition & 1 deletion Cassiopee/XCore/XCore/AdaptMesh/Pyra.h
Original file line number Diff line number Diff line change
Expand Up @@ -18,7 +18,7 @@
*/
#pragma once

#include "../common/common.h"
#include "common/common.h"

extern const Int normalIn_Py[5];

Expand Down
2 changes: 1 addition & 1 deletion Cassiopee/XCore/XCore/AdaptMesh/Tetra.h
Original file line number Diff line number Diff line change
Expand Up @@ -18,7 +18,7 @@
*/
#pragma once

#include "../common/common.h"
#include "common/common.h"

extern const Int normalIn_T[4];

Expand Down
3 changes: 1 addition & 2 deletions Cassiopee/XCore/XCore/AdaptMesh/Tri.h
Original file line number Diff line number Diff line change
Expand Up @@ -18,8 +18,7 @@
*/
#pragma once

#include "../common/common.h"

#include "common/common.h"
#include "Mesh.h"

void refine_tri(Int tri, Mesh *M);
Expand Down
2 changes: 1 addition & 1 deletion Cassiopee/XCore/XCore/intersectMesh/meshRefine.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -19,7 +19,7 @@
#include <stack>

#include "mesh.h"
#include "../common/common.h"
#include "common/common.h"

int meshes_mutual_refinement(IMesh &M, IMesh &S)
{
Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -19,7 +19,7 @@
#include <unordered_map>

#include "xcore.h"
#include "../common/common.h"
#include "common/common.h"
#include "karray.h"
#include "mesh.h"
#include "ray.h"
Expand Down
Loading

0 comments on commit 95eb4ec

Please sign in to comment.