2009-01-26 00:37:14 +05:00
|
|
|
#ifdef PARALLEL
|
|
|
|
|
|
|
|
#include <meshing.hpp>
|
|
|
|
#include "paralleltop.hpp"
|
|
|
|
|
2012-08-20 20:10:23 +06:00
|
|
|
// #define METIS4
|
2011-11-01 18:54:07 +06:00
|
|
|
|
2009-01-26 00:37:14 +05:00
|
|
|
|
|
|
|
#ifdef METIS
|
2011-07-15 03:36:19 +06:00
|
|
|
namespace metis {
|
2011-02-14 17:27:18 +05:00
|
|
|
extern "C" {
|
2011-11-01 18:54:07 +06:00
|
|
|
|
2012-11-01 18:57:40 +06:00
|
|
|
#include <metis.h>
|
2012-08-20 20:10:23 +06:00
|
|
|
|
|
|
|
#if METIS_VER_MAJOR >= 5
|
|
|
|
#define METIS5
|
2011-11-01 18:54:07 +06:00
|
|
|
typedef idx_t idxtype;
|
2012-08-20 20:10:23 +06:00
|
|
|
#else
|
|
|
|
#define METIS4
|
2012-09-28 18:39:41 +06:00
|
|
|
typedef idxtype idx_t;
|
2011-11-01 18:54:07 +06:00
|
|
|
#endif
|
2011-02-14 17:27:18 +05:00
|
|
|
}
|
|
|
|
}
|
2011-11-01 18:54:07 +06:00
|
|
|
|
2011-02-14 17:27:18 +05:00
|
|
|
using namespace metis;
|
|
|
|
#endif
|
2009-01-26 00:37:14 +05:00
|
|
|
|
|
|
|
namespace netgen
|
|
|
|
{
|
|
|
|
|
2011-06-30 18:49:38 +06:00
|
|
|
template <>
|
2012-08-20 20:10:23 +06:00
|
|
|
inline MPI_Datatype MyGetMPIType<PointIndex> ( )
|
2011-06-30 18:49:38 +06:00
|
|
|
{ return MPI_INT; }
|
|
|
|
|
|
|
|
|
2011-07-04 18:29:02 +06:00
|
|
|
void Mesh :: SendRecvMesh ()
|
|
|
|
{
|
2019-02-13 02:11:55 +05:00
|
|
|
int id = GetCommunicator().Rank();
|
|
|
|
int np = GetCommunicator().Size();
|
2011-07-04 18:29:02 +06:00
|
|
|
|
2018-07-19 20:33:26 +05:00
|
|
|
if (np == 1) {
|
|
|
|
throw NgException("SendRecvMesh called, but only one rank in communicator!!");
|
2011-07-04 18:29:02 +06:00
|
|
|
}
|
2018-07-19 20:33:26 +05:00
|
|
|
|
|
|
|
if (id == 0)
|
|
|
|
PrintMessage (1, "Send/Receive mesh");
|
2011-07-04 18:29:02 +06:00
|
|
|
|
2018-07-20 17:00:34 +05:00
|
|
|
// Why is this here??
|
2018-07-19 20:33:26 +05:00
|
|
|
if (id == 0)
|
|
|
|
{
|
|
|
|
paralleltop -> SetNV (GetNV());
|
|
|
|
paralleltop -> SetNE (GetNE());
|
|
|
|
paralleltop -> SetNSegm (GetNSeg());
|
|
|
|
paralleltop -> SetNSE (GetNSE());
|
|
|
|
}
|
2011-07-04 18:29:02 +06:00
|
|
|
|
|
|
|
if (id == 0)
|
|
|
|
SendMesh ();
|
|
|
|
else
|
|
|
|
ReceiveParallelMesh();
|
2012-06-21 20:19:25 +06:00
|
|
|
|
|
|
|
paralleltop -> UpdateCoarseGrid();
|
2011-07-15 03:36:19 +06:00
|
|
|
}
|
2011-07-14 00:32:11 +06:00
|
|
|
|
|
|
|
|
2011-07-04 18:29:02 +06:00
|
|
|
|
2009-01-26 00:37:14 +05:00
|
|
|
|
|
|
|
|
|
|
|
|
2011-06-30 18:49:38 +06:00
|
|
|
|
2011-07-15 03:36:19 +06:00
|
|
|
void Mesh :: SendMesh () const
|
|
|
|
{
|
|
|
|
Array<MPI_Request> sendrequests;
|
2019-01-30 15:27:21 +05:00
|
|
|
|
2019-02-13 02:11:55 +05:00
|
|
|
NgMPI_Comm comm = GetCommunicator();
|
|
|
|
int id = comm.Rank();
|
|
|
|
int ntasks = comm.Size();
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2018-07-19 20:33:26 +05:00
|
|
|
int dim = GetDimension();
|
2019-01-29 20:10:02 +05:00
|
|
|
MyMPI_Bcast(dim, comm);
|
2011-07-17 05:13:26 +06:00
|
|
|
|
2018-07-19 20:33:26 +05:00
|
|
|
|
2019-02-28 15:58:22 +05:00
|
|
|
// If the topology is not already updated, we do not need to
|
|
|
|
// build edges/faces.
|
|
|
|
auto & top = const_cast<MeshTopology&>(GetTopology());
|
|
|
|
if(top.NeedsUpdate()) {
|
|
|
|
top.SetBuildEdges(false);
|
|
|
|
top.SetBuildFaces(false);
|
|
|
|
top.Update();
|
|
|
|
}
|
2018-07-19 20:33:26 +05:00
|
|
|
|
|
|
|
PrintMessage ( 3, "Sending nr of elements");
|
|
|
|
|
2011-07-15 03:36:19 +06:00
|
|
|
Array<int> num_els_on_proc(ntasks);
|
|
|
|
num_els_on_proc = 0;
|
|
|
|
for (ElementIndex ei = 0; ei < GetNE(); ei++)
|
2019-02-15 18:34:27 +05:00
|
|
|
// num_els_on_proc[(*this)[ei].GetPartition()]++;
|
|
|
|
num_els_on_proc[vol_partition[ei]]++;
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2018-07-19 20:33:26 +05:00
|
|
|
MPI_Scatter (&num_els_on_proc[0], 1, MPI_INT,
|
|
|
|
MPI_IN_PLACE, -1, MPI_INT, 0, comm);
|
|
|
|
|
2011-07-15 03:36:19 +06:00
|
|
|
TABLE<ElementIndex> els_of_proc (num_els_on_proc);
|
|
|
|
for (ElementIndex ei = 0; ei < GetNE(); ei++)
|
2019-02-15 18:34:27 +05:00
|
|
|
// els_of_proc.Add ( (*this)[ei].GetPartition(), ei);
|
|
|
|
els_of_proc.Add (vol_partition[ei], ei);
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2018-07-19 20:33:26 +05:00
|
|
|
PrintMessage ( 3, "Building vertex/proc mapping");
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2011-07-17 05:13:26 +06:00
|
|
|
Array<int> num_sels_on_proc(ntasks);
|
|
|
|
num_sels_on_proc = 0;
|
|
|
|
for (SurfaceElementIndex ei = 0; ei < GetNSE(); ei++)
|
2019-02-15 18:34:27 +05:00
|
|
|
// num_sels_on_proc[(*this)[ei].GetPartition()]++;
|
|
|
|
num_sels_on_proc[surf_partition[ei]]++;
|
2011-07-17 05:13:26 +06:00
|
|
|
|
|
|
|
TABLE<SurfaceElementIndex> sels_of_proc (num_sels_on_proc);
|
|
|
|
for (SurfaceElementIndex ei = 0; ei < GetNSE(); ei++)
|
2019-02-15 18:34:27 +05:00
|
|
|
// sels_of_proc.Add ( (*this)[ei].GetPartition(), ei);
|
|
|
|
sels_of_proc.Add (surf_partition[ei], ei);
|
2011-07-17 05:13:26 +06:00
|
|
|
|
2012-08-20 20:10:23 +06:00
|
|
|
Array<int> num_segs_on_proc(ntasks);
|
|
|
|
num_segs_on_proc = 0;
|
|
|
|
for (SegmentIndex ei = 0; ei < GetNSeg(); ei++)
|
2019-02-15 18:34:27 +05:00
|
|
|
// num_segs_on_proc[(*this)[ei].GetPartition()]++;
|
|
|
|
num_segs_on_proc[seg_partition[ei]]++;
|
2012-08-20 20:10:23 +06:00
|
|
|
|
|
|
|
TABLE<SegmentIndex> segs_of_proc (num_segs_on_proc);
|
|
|
|
for (SegmentIndex ei = 0; ei < GetNSeg(); ei++)
|
2019-02-15 18:34:27 +05:00
|
|
|
segs_of_proc.Add (seg_partition[ei], ei);
|
2012-08-20 20:10:23 +06:00
|
|
|
|
2011-07-17 05:13:26 +06:00
|
|
|
|
2018-07-19 20:33:26 +05:00
|
|
|
/**
|
2018-07-20 16:44:46 +05:00
|
|
|
----- STRATEGY FOR PERIODIC MESHES -----
|
|
|
|
|
2018-07-19 20:33:26 +05:00
|
|
|
Whenever two vertices are identified by periodicity, any proc
|
|
|
|
that gets one of the vertices actually gets both of them.
|
|
|
|
This has to be transitive, that is, if
|
|
|
|
a <-> b and b <-> c,
|
|
|
|
then any proc that has vertex a also has vertices b and c!
|
2018-07-20 16:44:46 +05:00
|
|
|
|
|
|
|
Surfaceelements and Segments that are identified by
|
|
|
|
periodicity are treated the same way.
|
|
|
|
|
|
|
|
We need to duplicate these so we have containers to
|
|
|
|
hold the edges/facets. Afaik, a mesh cannot have nodes
|
|
|
|
that are not part of some sort of element.
|
|
|
|
|
2018-07-19 20:33:26 +05:00
|
|
|
**/
|
2011-07-17 05:13:26 +06:00
|
|
|
|
2018-07-20 16:44:46 +05:00
|
|
|
/** First, we build tables for vertex identification. **/
|
2018-07-19 20:33:26 +05:00
|
|
|
Array<INDEX_2> per_pairs;
|
|
|
|
Array<INDEX_2> pp2;
|
|
|
|
auto & idents = GetIdentifications();
|
|
|
|
bool has_periodic = false;
|
|
|
|
for (int idnr = 1; idnr < idents.GetMaxNr()+1; idnr++)
|
2011-07-15 03:36:19 +06:00
|
|
|
{
|
2018-07-19 20:33:26 +05:00
|
|
|
if(idents.GetType(idnr)!=Identifications::PERIODIC) continue;
|
|
|
|
has_periodic = true;
|
|
|
|
idents.GetPairs(idnr, pp2);
|
|
|
|
per_pairs.Append(pp2);
|
|
|
|
}
|
|
|
|
Array<int, PointIndex::BASE> npvs(GetNV());
|
|
|
|
npvs = 0;
|
|
|
|
for (int k = 0; k < per_pairs.Size(); k++) {
|
|
|
|
npvs[per_pairs[k].I1()]++;
|
|
|
|
npvs[per_pairs[k].I2()]++;
|
|
|
|
}
|
2018-07-20 16:44:46 +05:00
|
|
|
|
2018-07-19 20:33:26 +05:00
|
|
|
/** for each vertex, gives us all identified vertices **/
|
|
|
|
TABLE<PointIndex, PointIndex::BASE> per_verts(npvs);
|
|
|
|
for (int k = 0; k < per_pairs.Size(); k++) {
|
|
|
|
per_verts.Add(per_pairs[k].I1(), per_pairs[k].I2());
|
|
|
|
per_verts.Add(per_pairs[k].I2(), per_pairs[k].I1());
|
|
|
|
}
|
|
|
|
for (int k = PointIndex::BASE; k < GetNV()+PointIndex::BASE; k++) {
|
|
|
|
BubbleSort(per_verts[k]);
|
|
|
|
}
|
2018-07-20 16:44:46 +05:00
|
|
|
|
|
|
|
/** The same table as per_verts, but TRANSITIVE!! **/
|
2018-07-19 20:33:26 +05:00
|
|
|
auto iterate_per_verts_trans = [&](auto f){
|
|
|
|
Array<int> allvs;
|
|
|
|
for (int k = PointIndex::BASE; k < GetNV()+PointIndex::BASE; k++)
|
|
|
|
{
|
|
|
|
allvs.SetSize(0);
|
|
|
|
allvs.Append(per_verts[k]);
|
|
|
|
bool changed = true;
|
|
|
|
while(changed) {
|
|
|
|
changed = false;
|
|
|
|
for (int j = 0; j<allvs.Size(); j++)
|
2011-07-17 05:13:26 +06:00
|
|
|
{
|
2018-07-19 20:33:26 +05:00
|
|
|
auto pervs2 = per_verts[allvs[j]];
|
|
|
|
for (int l = 0; l < pervs2.Size(); l++)
|
2011-07-17 05:13:26 +06:00
|
|
|
{
|
2018-07-19 20:33:26 +05:00
|
|
|
auto addv = pervs2[l];
|
|
|
|
if (allvs.Contains(addv) || addv==k) continue;
|
|
|
|
changed = true;
|
|
|
|
allvs.Append(addv);
|
2012-08-20 20:10:23 +06:00
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
2018-07-19 20:33:26 +05:00
|
|
|
f(k, allvs);
|
|
|
|
}
|
|
|
|
};
|
|
|
|
iterate_per_verts_trans([&](auto k, auto & allvs) {
|
|
|
|
npvs[k] = allvs.Size();
|
|
|
|
});
|
|
|
|
TABLE<PointIndex, PointIndex::BASE> per_verts_trans(npvs);
|
|
|
|
iterate_per_verts_trans([&](auto k, auto & allvs) {
|
|
|
|
for (int j = 0; j<allvs.Size(); j++)
|
|
|
|
per_verts_trans.Add(k, allvs[j]);
|
|
|
|
});
|
|
|
|
for (int k = PointIndex::BASE; k < GetNV()+PointIndex::BASE; k++) {
|
|
|
|
BubbleSort(per_verts_trans[k]);
|
|
|
|
}
|
2012-08-20 20:10:23 +06:00
|
|
|
|
2018-07-20 16:44:46 +05:00
|
|
|
/** Now we build the vertex-data to send to the workers. **/
|
|
|
|
Array<int, PointIndex::BASE> vert_flag (GetNV());
|
|
|
|
Array<int, PointIndex::BASE> num_procs_on_vert (GetNV());
|
|
|
|
Array<int> num_verts_on_proc (ntasks);
|
|
|
|
num_verts_on_proc = 0;
|
|
|
|
num_procs_on_vert = 0;
|
2018-07-19 20:33:26 +05:00
|
|
|
auto iterate_vertices = [&](auto f) {
|
|
|
|
vert_flag = -1;
|
|
|
|
for (int dest = 1; dest < ntasks; dest++)
|
|
|
|
{
|
|
|
|
FlatArray<ElementIndex> els = els_of_proc[dest];
|
|
|
|
for (int hi = 0; hi < els.Size(); hi++)
|
|
|
|
{
|
|
|
|
const Element & el = (*this) [ els[hi] ];
|
|
|
|
for (int i = 0; i < el.GetNP(); i++)
|
|
|
|
f(el[i], dest);
|
|
|
|
}
|
|
|
|
FlatArray<SurfaceElementIndex> sels = sels_of_proc[dest];
|
|
|
|
for (int hi = 0; hi < sels.Size(); hi++)
|
|
|
|
{
|
|
|
|
const Element2d & el = (*this) [ sels[hi] ];
|
|
|
|
for (int i = 0; i < el.GetNP(); i++)
|
|
|
|
f(el[i], dest);
|
|
|
|
}
|
|
|
|
FlatArray<SegmentIndex> segs = segs_of_proc[dest];
|
|
|
|
for (int hi = 0; hi < segs.Size(); hi++)
|
|
|
|
{
|
|
|
|
const Segment & el = (*this) [segs[hi]];
|
|
|
|
for (int i = 0; i < 2; i++)
|
|
|
|
f(el[i], dest);
|
|
|
|
}
|
|
|
|
}
|
|
|
|
};
|
|
|
|
/** count vertices per proc and procs per vertex **/
|
|
|
|
iterate_vertices([&](auto vertex, auto dest){
|
|
|
|
auto countit = [&] (auto vertex, auto dest) {
|
|
|
|
if (vert_flag[vertex] < dest)
|
|
|
|
{
|
|
|
|
vert_flag[vertex] = dest;
|
|
|
|
num_verts_on_proc[dest]++;
|
|
|
|
num_procs_on_vert[vertex]++;
|
|
|
|
GetParallelTopology().SetDistantPNum (dest, vertex);
|
|
|
|
}
|
|
|
|
};
|
|
|
|
countit(vertex, dest);
|
|
|
|
auto pers = per_verts_trans[vertex];
|
|
|
|
for(int j = 0; j < pers.Size(); j++)
|
|
|
|
countit(pers[j], dest);
|
|
|
|
});
|
2011-07-15 03:36:19 +06:00
|
|
|
TABLE<PointIndex> verts_of_proc (num_verts_on_proc);
|
|
|
|
TABLE<int, PointIndex::BASE> procs_of_vert (num_procs_on_vert);
|
|
|
|
TABLE<int, PointIndex::BASE> loc_num_of_vert (num_procs_on_vert);
|
2018-07-19 20:33:26 +05:00
|
|
|
/** Write vertex/proc mappingfs to tables **/
|
|
|
|
iterate_vertices([&](auto vertex, auto dest) {
|
|
|
|
auto addit = [&] (auto vertex, auto dest) {
|
|
|
|
if (vert_flag[vertex] < dest)
|
|
|
|
{
|
|
|
|
vert_flag[vertex] = dest;
|
|
|
|
procs_of_vert.Add (vertex, dest);
|
|
|
|
}
|
|
|
|
};
|
|
|
|
addit(vertex, dest);
|
|
|
|
auto pers = per_verts_trans[vertex];
|
|
|
|
for(int j = 0; j < pers.Size(); j++)
|
|
|
|
addit(pers[j], dest);
|
|
|
|
});
|
2018-07-20 16:44:46 +05:00
|
|
|
/**
|
|
|
|
local vertex numbers on distant procs
|
|
|
|
(I think this was only used for debugging??)
|
|
|
|
**/
|
2011-07-15 03:36:19 +06:00
|
|
|
for (int vert = 1; vert <= GetNP(); vert++ )
|
|
|
|
{
|
|
|
|
FlatArray<int> procs = procs_of_vert[vert];
|
|
|
|
for (int j = 0; j < procs.Size(); j++)
|
|
|
|
{
|
|
|
|
int dest = procs[j];
|
|
|
|
verts_of_proc.Add (dest, vert);
|
|
|
|
loc_num_of_vert.Add (vert, verts_of_proc[dest].Size());
|
|
|
|
}
|
|
|
|
}
|
2018-07-19 20:33:26 +05:00
|
|
|
PrintMessage ( 3, "Sending Vertices - vertices");
|
|
|
|
|
2011-07-15 03:36:19 +06:00
|
|
|
for (int dest = 1; dest < ntasks; dest++)
|
|
|
|
{
|
|
|
|
FlatArray<PointIndex> verts = verts_of_proc[dest];
|
2019-01-29 20:10:02 +05:00
|
|
|
sendrequests.Append (MyMPI_ISend (verts, dest, MPI_TAG_MESH+1, comm));
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2011-07-15 03:36:19 +06:00
|
|
|
MPI_Datatype mptype = MeshPoint::MyGetMPIType();
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2011-07-15 03:36:19 +06:00
|
|
|
int numv = verts.Size();
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2011-07-15 03:36:19 +06:00
|
|
|
MPI_Datatype newtype;
|
2012-06-13 15:07:11 +06:00
|
|
|
Array<int> blocklen (numv);
|
2011-07-15 03:36:19 +06:00
|
|
|
blocklen = 1;
|
|
|
|
|
|
|
|
MPI_Type_indexed (numv, &blocklen[0],
|
|
|
|
reinterpret_cast<int*> (&verts[0]),
|
|
|
|
mptype, &newtype);
|
|
|
|
MPI_Type_commit (&newtype);
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2011-07-15 03:36:19 +06:00
|
|
|
MPI_Request request;
|
2019-01-29 20:10:02 +05:00
|
|
|
MPI_Isend( &points[0], 1, newtype, dest, MPI_TAG_MESH+1, comm, &request);
|
2011-07-15 03:36:19 +06:00
|
|
|
sendrequests.Append (request);
|
|
|
|
}
|
2009-01-26 00:37:14 +05:00
|
|
|
|
|
|
|
Array<int> num_distpnums(ntasks);
|
|
|
|
num_distpnums = 0;
|
|
|
|
|
2018-07-19 20:33:26 +05:00
|
|
|
|
|
|
|
/**
|
2018-07-20 16:44:46 +05:00
|
|
|
Next, we send the identifications themselfs.
|
|
|
|
|
2018-07-19 20:33:26 +05:00
|
|
|
Info about periodic identifications sent to each proc is an array of
|
|
|
|
integers.
|
|
|
|
- maxidentnr
|
|
|
|
- type for each identification
|
|
|
|
- nr of pairs for each identification (each pair is local!)
|
|
|
|
- pairs for each periodic ident (global numbers)
|
|
|
|
**/
|
2018-07-20 16:44:46 +05:00
|
|
|
PrintMessage ( 3, "Sending Vertices - identifications");
|
2018-07-19 20:33:26 +05:00
|
|
|
int maxidentnr = idents.GetMaxNr();
|
|
|
|
Array<int> ppd_sizes(ntasks);
|
|
|
|
ppd_sizes = 1 + 2*maxidentnr;
|
|
|
|
for (int idnr = 1; idnr < idents.GetMaxNr()+1; idnr++)
|
|
|
|
{
|
|
|
|
if(idents.GetType(idnr)!=Identifications::PERIODIC) continue;
|
|
|
|
idents.GetPairs(idnr, pp2);
|
|
|
|
for(int j = 0; j<pp2.Size(); j++)
|
|
|
|
{
|
|
|
|
INDEX_2 & pair = pp2[j];
|
|
|
|
// both are on same procs!
|
|
|
|
auto ps = procs_of_vert[pair.I1()];
|
|
|
|
for (int l = 0; l < ps.Size(); l++)
|
|
|
|
{
|
|
|
|
ppd_sizes[ps[l]] += 2;
|
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
|
|
|
TABLE<int> pp_data(ppd_sizes);
|
|
|
|
for(int dest = 0; dest < ntasks; dest++)
|
|
|
|
pp_data.Add(dest, maxidentnr);
|
|
|
|
for (int dest = 0; dest < ntasks; dest++)
|
|
|
|
{
|
|
|
|
for (int idnr = 1; idnr < idents.GetMaxNr()+1; idnr++)
|
|
|
|
pp_data.Add(dest, idents.GetType(idnr));
|
|
|
|
for (int idnr = 1; idnr < idents.GetMaxNr()+1; idnr++)
|
|
|
|
pp_data.Add(dest, 0);
|
|
|
|
}
|
|
|
|
for (int idnr = 1; idnr < idents.GetMaxNr()+1; idnr++)
|
|
|
|
{
|
|
|
|
if(idents.GetType(idnr)!=Identifications::PERIODIC) continue;
|
|
|
|
idents.GetPairs(idnr, pp2);
|
|
|
|
for(int j = 0; j<pp2.Size(); j++)
|
|
|
|
{
|
|
|
|
INDEX_2 & pair = pp2[j];
|
|
|
|
auto ps = procs_of_vert[pair.I1()];
|
|
|
|
for (int l = 0; l < ps.Size(); l++)
|
|
|
|
{
|
|
|
|
auto p = ps[l];
|
|
|
|
pp_data[p][maxidentnr + idnr]++;
|
|
|
|
pp_data.Add(p, pair.I1());
|
|
|
|
pp_data.Add(p, pair.I2());
|
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
|
|
|
Array<MPI_Request> req_per;
|
|
|
|
for(int dest = 1; dest < ntasks; dest++)
|
2019-01-29 20:10:02 +05:00
|
|
|
req_per.Append(MyMPI_ISend(pp_data[dest], dest, MPI_TAG_MESH+1, comm));
|
2018-07-19 20:33:26 +05:00
|
|
|
MPI_Waitall(req_per.Size(), &req_per[0], MPI_STATUS_IGNORE);
|
2018-07-20 16:44:46 +05:00
|
|
|
|
2018-07-19 20:33:26 +05:00
|
|
|
PrintMessage ( 3, "Sending Vertices - distprocs");
|
|
|
|
|
2011-07-14 00:32:11 +06:00
|
|
|
for (int vert = 1; vert <= GetNP(); vert++)
|
2009-01-26 00:37:14 +05:00
|
|
|
{
|
|
|
|
FlatArray<int> procs = procs_of_vert[vert];
|
|
|
|
for (int j = 0; j < procs.Size(); j++)
|
|
|
|
num_distpnums[procs[j]] += 3 * (procs.Size()-1);
|
|
|
|
}
|
|
|
|
|
|
|
|
TABLE<int> distpnums (num_distpnums);
|
|
|
|
|
2011-07-14 00:32:11 +06:00
|
|
|
for (int vert = 1; vert <= GetNP(); vert++)
|
2009-01-26 00:37:14 +05:00
|
|
|
{
|
|
|
|
FlatArray<int> procs = procs_of_vert[vert];
|
|
|
|
for (int j = 0; j < procs.Size(); j++)
|
|
|
|
for (int k = 0; k < procs.Size(); k++)
|
|
|
|
if (j != k)
|
|
|
|
{
|
|
|
|
distpnums.Add (procs[j], loc_num_of_vert[vert][j]);
|
|
|
|
distpnums.Add (procs[j], procs_of_vert[vert][k]);
|
|
|
|
distpnums.Add (procs[j], loc_num_of_vert[vert][k]);
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
for ( int dest = 1; dest < ntasks; dest ++ )
|
2019-01-29 20:10:02 +05:00
|
|
|
sendrequests.Append (MyMPI_ISend (distpnums[dest], dest, MPI_TAG_MESH+1, comm));
|
2011-06-30 18:49:38 +06:00
|
|
|
|
2011-07-17 05:13:26 +06:00
|
|
|
|
|
|
|
|
2009-01-26 00:37:14 +05:00
|
|
|
PrintMessage ( 3, "Sending elements" );
|
|
|
|
|
2011-07-04 18:29:02 +06:00
|
|
|
Array<int> elarraysize (ntasks);
|
|
|
|
elarraysize = 0;
|
2011-07-14 00:32:11 +06:00
|
|
|
for ( int ei = 1; ei <= GetNE(); ei++)
|
2011-07-04 18:29:02 +06:00
|
|
|
{
|
2011-07-14 00:32:11 +06:00
|
|
|
const Element & el = VolumeElement (ei);
|
2019-02-15 18:34:27 +05:00
|
|
|
// int dest = el.GetPartition();
|
|
|
|
int dest = vol_partition[ei-1];
|
2011-07-04 18:29:02 +06:00
|
|
|
elarraysize[dest] += 3 + el.GetNP();
|
|
|
|
}
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2011-07-04 18:29:02 +06:00
|
|
|
TABLE<int> elementarrays(elarraysize);
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2011-07-14 00:32:11 +06:00
|
|
|
for (int ei = 1; ei <= GetNE(); ei++)
|
2009-01-26 00:37:14 +05:00
|
|
|
{
|
2011-07-14 00:32:11 +06:00
|
|
|
const Element & el = VolumeElement (ei);
|
2019-02-15 18:34:27 +05:00
|
|
|
// int dest = el.GetPartition();
|
|
|
|
int dest = vol_partition[ei-1];
|
|
|
|
|
2011-07-04 18:29:02 +06:00
|
|
|
elementarrays.Add (dest, ei);
|
|
|
|
elementarrays.Add (dest, el.GetIndex());
|
|
|
|
elementarrays.Add (dest, el.GetNP());
|
|
|
|
for (int i = 0; i < el.GetNP(); i++)
|
|
|
|
elementarrays.Add (dest, el[i]);
|
2009-01-26 00:37:14 +05:00
|
|
|
}
|
|
|
|
|
2011-07-04 18:29:02 +06:00
|
|
|
for (int dest = 1; dest < ntasks; dest ++ )
|
2019-01-29 20:10:02 +05:00
|
|
|
sendrequests.Append (MyMPI_ISend (elementarrays[dest], dest, MPI_TAG_MESH+2, comm));
|
2011-07-04 18:29:02 +06:00
|
|
|
|
2011-08-22 20:51:21 +06:00
|
|
|
|
2009-01-26 00:37:14 +05:00
|
|
|
PrintMessage ( 3, "Sending Face Descriptors" );
|
|
|
|
|
2011-07-14 00:32:11 +06:00
|
|
|
Array<double> fddata (6 * GetNFD());
|
|
|
|
for (int fdi = 1; fdi <= GetNFD(); fdi++)
|
|
|
|
{
|
|
|
|
fddata[6*fdi-6] = GetFaceDescriptor(fdi).SurfNr();
|
|
|
|
fddata[6*fdi-5] = GetFaceDescriptor(fdi).DomainIn();
|
|
|
|
fddata[6*fdi-4] = GetFaceDescriptor(fdi).DomainOut();
|
|
|
|
fddata[6*fdi-3] = GetFaceDescriptor(fdi).BCProperty();
|
|
|
|
fddata[6*fdi-2] = GetFaceDescriptor(fdi).domin_singular;
|
|
|
|
fddata[6*fdi-1] = GetFaceDescriptor(fdi).domout_singular;
|
|
|
|
|
|
|
|
}
|
|
|
|
for (int dest = 1; dest < ntasks; dest++)
|
2019-01-29 20:10:02 +05:00
|
|
|
sendrequests.Append (MyMPI_ISend (fddata, dest, MPI_TAG_MESH+3, comm));
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2018-07-19 20:33:26 +05:00
|
|
|
/** Surface Elements **/
|
|
|
|
|
2009-01-26 00:37:14 +05:00
|
|
|
PrintMessage ( 3, "Sending Surface elements" );
|
2018-07-19 20:33:26 +05:00
|
|
|
// build sel-identification
|
|
|
|
size_t nse = GetNSE();
|
|
|
|
Array<SurfaceElementIndex> ided_sel(nse);
|
|
|
|
ided_sel = -1;
|
|
|
|
bool has_ided_sels = false;
|
|
|
|
if(GetNE() && has_periodic) //we can only have identified surf-els if we have vol-els (right?)
|
|
|
|
{
|
|
|
|
Array<SurfaceElementIndex, 0> os1, os2;
|
|
|
|
for(SurfaceElementIndex sei = 0; sei < GetNSE(); sei++)
|
|
|
|
{
|
|
|
|
if(ided_sel[sei]!=-1) continue;
|
|
|
|
const Element2d & sel = (*this)[sei];
|
|
|
|
FlatArray<const PointIndex> points = sel.PNums();
|
|
|
|
auto ided1 = per_verts[points[0]];
|
|
|
|
os1.SetSize(0);
|
|
|
|
for (int j = 0; j < ided1.Size(); j++)
|
|
|
|
os1.Append(GetTopology().GetVertexSurfaceElements(ided1[j]));
|
|
|
|
for (int j = 1; j < points.Size(); j++)
|
|
|
|
{
|
|
|
|
os2.SetSize(0);
|
|
|
|
auto p2 = points[j];
|
|
|
|
auto ided2 = per_verts[p2];
|
|
|
|
for (int l = 0; l < ided2.Size(); l++)
|
|
|
|
os2.Append(GetTopology().GetVertexSurfaceElements(ided2[l]));
|
|
|
|
for (int m = 0; m<os1.Size(); m++) {
|
|
|
|
if(!os2.Contains(os1[m])) {
|
|
|
|
os1.Delete(m);
|
|
|
|
m--;
|
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
|
|
|
if(!os1.Size()) continue;
|
|
|
|
if(os1.Size()>1) {
|
|
|
|
throw NgException("SurfaceElement identified with more than one other??");
|
|
|
|
}
|
|
|
|
const Element2d & sel2 = (*this)[sei];
|
|
|
|
FlatArray<const PointIndex> points2 = sel2.PNums();
|
|
|
|
has_ided_sels = true;
|
|
|
|
ided_sel[sei] = os1[0];
|
|
|
|
ided_sel[os1[0]] = sei;
|
|
|
|
}
|
|
|
|
}
|
|
|
|
// build sel data to send
|
|
|
|
auto iterate_sels = [&](auto f) {
|
|
|
|
for (SurfaceElementIndex sei = 0; sei < GetNSE(); sei++ )
|
|
|
|
{
|
|
|
|
const Element2d & sel = (*this)[sei];
|
2019-02-15 18:34:27 +05:00
|
|
|
// int dest = (*this)[sei].GetPartition();
|
|
|
|
int dest = surf_partition[sei];
|
2018-07-19 20:33:26 +05:00
|
|
|
f(sei, sel, dest);
|
|
|
|
if(ided_sel[sei]!=-1)
|
|
|
|
{
|
2019-02-15 18:34:27 +05:00
|
|
|
// int dest2 = (*this)[ided_sel[sei]].GetPartition();
|
|
|
|
int dest2 = surf_partition[ided_sel[sei]];
|
2018-07-19 20:33:26 +05:00
|
|
|
f(sei, sel, dest2);
|
|
|
|
}
|
|
|
|
}
|
|
|
|
};
|
|
|
|
Array <int> nlocsel(ntasks), bufsize(ntasks);
|
2011-07-17 05:13:26 +06:00
|
|
|
nlocsel = 0;
|
|
|
|
bufsize = 1;
|
2018-07-19 20:33:26 +05:00
|
|
|
iterate_sels([&](SurfaceElementIndex sei, const Element2d & sel, int dest){
|
|
|
|
nlocsel[dest]++;
|
2011-07-17 05:13:26 +06:00
|
|
|
bufsize[dest] += 4 + 2*sel.GetNP();
|
2018-07-19 20:33:26 +05:00
|
|
|
});
|
2011-07-14 00:32:11 +06:00
|
|
|
TABLE<int> selbuf(bufsize);
|
|
|
|
for (int dest = 1; dest < ntasks; dest++ )
|
|
|
|
selbuf.Add (dest, nlocsel[dest]);
|
2018-07-19 20:33:26 +05:00
|
|
|
iterate_sels([&](SurfaceElementIndex sei, const auto & sel, int dest) {
|
2011-07-17 05:13:26 +06:00
|
|
|
selbuf.Add (dest, sei);
|
|
|
|
selbuf.Add (dest, sel.GetIndex());
|
2014-04-28 13:09:39 +06:00
|
|
|
// selbuf.Add (dest, 0);
|
2011-07-17 05:13:26 +06:00
|
|
|
selbuf.Add (dest, sel.GetNP());
|
|
|
|
for ( int ii = 1; ii <= sel.GetNP(); ii++)
|
|
|
|
{
|
|
|
|
selbuf.Add (dest, sel.PNum(ii));
|
|
|
|
selbuf.Add (dest, sel.GeomInfoPi(ii).trignum);
|
2009-01-26 00:37:14 +05:00
|
|
|
}
|
2018-07-19 20:33:26 +05:00
|
|
|
});
|
|
|
|
// distribute sel data
|
2011-07-14 00:32:11 +06:00
|
|
|
for (int dest = 1; dest < ntasks; dest++)
|
2019-01-29 20:10:02 +05:00
|
|
|
sendrequests.Append (MyMPI_ISend(selbuf[dest], dest, MPI_TAG_MESH+4, comm));
|
2018-07-19 20:33:26 +05:00
|
|
|
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2018-07-20 16:44:46 +05:00
|
|
|
/** Segments **/
|
2009-01-26 00:37:14 +05:00
|
|
|
PrintMessage ( 3, "Sending Edge Segments");
|
2018-07-19 20:33:26 +05:00
|
|
|
auto iterate_segs1 = [&](auto f) {
|
|
|
|
Array<SegmentIndex> osegs1, osegs2, osegs_both;
|
|
|
|
Array<int> type1, type2;
|
|
|
|
for(SegmentIndex segi = 0; segi < GetNSeg(); segi++)
|
|
|
|
{
|
|
|
|
const Segment & seg = (*this)[segi];
|
|
|
|
int segnp = seg.GetNP();
|
|
|
|
PointIndex pi1 = seg[0];
|
|
|
|
auto ided1 = per_verts[pi1];
|
|
|
|
PointIndex pi2 = seg[1];
|
|
|
|
auto ided2 = per_verts[pi2];
|
|
|
|
if (!(ided1.Size() && ided2.Size())) continue;
|
|
|
|
osegs1.SetSize(0);
|
|
|
|
type1.SetSize(0);
|
|
|
|
for (int l = 0; l<ided1.Size(); l++)
|
|
|
|
{
|
|
|
|
auto ospart = GetTopology().GetVertexSegments(ided1[l]);
|
|
|
|
for(int j=0; j<ospart.Size(); j++)
|
|
|
|
{
|
|
|
|
if(osegs1.Contains(ospart[j]))
|
2018-07-20 16:44:46 +05:00
|
|
|
throw NgException("Periodic Mesh did something weird.");
|
2018-07-19 20:33:26 +05:00
|
|
|
osegs1.Append(ospart[j]);
|
|
|
|
type1.Append(idents.GetSymmetric(pi1, ided1[l]));
|
|
|
|
}
|
|
|
|
}
|
|
|
|
osegs2.SetSize(0);
|
|
|
|
type2.SetSize(0);
|
|
|
|
for (int l = 0; l<ided2.Size(); l++)
|
|
|
|
{
|
|
|
|
auto ospart = GetTopology().GetVertexSegments(ided2[l]);
|
|
|
|
for(int j=0; j<ospart.Size(); j++)
|
|
|
|
{
|
|
|
|
if(osegs2.Contains(ospart[j]))
|
2018-07-20 16:44:46 +05:00
|
|
|
throw NgException("Periodic Mesh did something weird.");
|
2018-07-19 20:33:26 +05:00
|
|
|
osegs2.Append(ospart[j]);
|
|
|
|
type2.Append(idents.GetSymmetric(pi2, ided2[l]));
|
|
|
|
}
|
|
|
|
}
|
|
|
|
osegs_both.SetSize(0);
|
|
|
|
for (int l = 0; l<osegs1.Size(); l++) {
|
|
|
|
auto pos = osegs2.Pos(osegs1[l]);
|
|
|
|
if (pos == -1) continue;
|
|
|
|
if (type1[l] != type2[pos]) continue;
|
|
|
|
osegs_both.Append(osegs1[l]);
|
|
|
|
}
|
|
|
|
for(int l = 0; l<osegs_both.Size(); l++) {
|
2019-02-04 17:12:31 +05:00
|
|
|
int segnp2 = (*this)[osegs_both[l]].GetNP();
|
|
|
|
if(segnp!=segnp2)
|
2018-07-19 20:33:26 +05:00
|
|
|
throw NgException("Tried to identify non-curved and curved Segment!");
|
|
|
|
}
|
|
|
|
for(int l = 0; l<osegs_both.Size(); l++) {
|
|
|
|
f(segi, osegs_both[l]);
|
|
|
|
}
|
|
|
|
}
|
|
|
|
};
|
|
|
|
Array<int> per_seg_size(GetNSeg());
|
|
|
|
per_seg_size = 0;
|
|
|
|
iterate_segs1([&](SegmentIndex segi1, SegmentIndex segi2)
|
|
|
|
{ per_seg_size[segi1]++; });
|
|
|
|
TABLE<SegmentIndex> per_seg(per_seg_size);
|
|
|
|
iterate_segs1([&](SegmentIndex segi1, SegmentIndex segi2)
|
|
|
|
{ per_seg.Add(segi1, segi2); });
|
|
|
|
// make per_seg transitive
|
|
|
|
auto iterate_per_seg_trans = [&](auto f){
|
|
|
|
Array<SegmentIndex> allsegs;
|
|
|
|
for (SegmentIndex segi = 0; segi < GetNSeg(); segi++)
|
|
|
|
{
|
|
|
|
allsegs.SetSize(0);
|
|
|
|
allsegs.Append(per_seg[segi]);
|
|
|
|
bool changed = true;
|
|
|
|
while (changed)
|
|
|
|
{
|
|
|
|
changed = false;
|
|
|
|
for (int j = 0; j<allsegs.Size(); j++)
|
|
|
|
{
|
|
|
|
auto persegs2 = per_seg[allsegs[j]];
|
|
|
|
for (int l = 0; l<persegs2.Size(); l++)
|
|
|
|
{
|
|
|
|
auto addseg = persegs2[l];
|
|
|
|
if (allsegs.Contains(addseg) || addseg==segi) continue;
|
|
|
|
allsegs.Append(addseg);
|
|
|
|
changed = true;
|
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
|
|
|
f(segi, allsegs);
|
|
|
|
}
|
|
|
|
};
|
|
|
|
iterate_per_seg_trans([&](SegmentIndex segi, Array<SegmentIndex> & segs){
|
|
|
|
for (int j = 0; j < segs.Size(); j++)
|
|
|
|
per_seg_size[segi] = segs.Size();
|
|
|
|
});
|
|
|
|
TABLE<SegmentIndex> per_seg_trans(per_seg_size);
|
|
|
|
iterate_per_seg_trans([&](SegmentIndex segi, Array<SegmentIndex> & segs){
|
|
|
|
for (int j = 0; j < segs.Size(); j++)
|
|
|
|
per_seg_trans.Add(segi, segs[j]);
|
|
|
|
});
|
|
|
|
// build segment data
|
|
|
|
Array<int> dests;
|
|
|
|
auto iterate_segs2 = [&](auto f)
|
|
|
|
{
|
|
|
|
for (SegmentIndex segi = 0; segi<GetNSeg(); segi++)
|
|
|
|
{
|
|
|
|
const Segment & seg = (*this)[segi];
|
|
|
|
dests.SetSize(0);
|
2019-02-15 18:34:27 +05:00
|
|
|
// dests.Append(seg.GetPartition());
|
|
|
|
dests.Append(seg_partition[segi]);
|
2018-07-19 20:33:26 +05:00
|
|
|
for (int l = 0; l < per_seg_trans[segi].Size(); l++)
|
|
|
|
{
|
2019-02-15 18:34:27 +05:00
|
|
|
// int dest2 = (*this)[per_seg_trans[segi][l]].GetPartition();
|
|
|
|
int dest2 = seg_partition[per_seg_trans[segi][l]];
|
2018-07-19 20:33:26 +05:00
|
|
|
if(!dests.Contains(dest2))
|
|
|
|
dests.Append(dest2);
|
|
|
|
}
|
|
|
|
for (int l = 0; l < dests.Size(); l++)
|
|
|
|
f(segi, seg, dests[l]);
|
|
|
|
}
|
|
|
|
};
|
|
|
|
Array<int> nloc_seg(ntasks);
|
|
|
|
// bufsize = 1; //was originally this - why??
|
|
|
|
bufsize = 0;
|
|
|
|
nloc_seg = 0;
|
|
|
|
iterate_segs2([&](auto segi, const auto & seg, int dest)
|
|
|
|
{
|
|
|
|
nloc_seg[dest]++;
|
|
|
|
bufsize[dest] += 14;
|
|
|
|
});
|
|
|
|
TABLE<double> segm_buf(bufsize);
|
|
|
|
iterate_segs2([&](auto segi, const auto & seg, int dest)
|
|
|
|
{
|
|
|
|
segm_buf.Add (dest, segi);
|
|
|
|
segm_buf.Add (dest, seg.si);
|
|
|
|
segm_buf.Add (dest, seg.pnums[0]);
|
|
|
|
segm_buf.Add (dest, seg.pnums[1]);
|
|
|
|
segm_buf.Add (dest, seg.geominfo[0].trignum);
|
|
|
|
segm_buf.Add (dest, seg.geominfo[1].trignum);
|
|
|
|
segm_buf.Add (dest, seg.surfnr1);
|
|
|
|
segm_buf.Add (dest, seg.surfnr2);
|
|
|
|
segm_buf.Add (dest, seg.edgenr);
|
|
|
|
segm_buf.Add (dest, seg.epgeominfo[0].dist);
|
|
|
|
segm_buf.Add (dest, seg.epgeominfo[1].edgenr);
|
|
|
|
segm_buf.Add (dest, seg.epgeominfo[1].dist);
|
|
|
|
segm_buf.Add (dest, seg.singedge_right);
|
|
|
|
segm_buf.Add (dest, seg.singedge_left);
|
|
|
|
});
|
|
|
|
// distrubute segment data
|
|
|
|
for (int dest = 1; dest < ntasks; dest++)
|
2019-01-29 20:10:02 +05:00
|
|
|
sendrequests.Append (MyMPI_ISend(segm_buf[dest], dest, MPI_TAG_MESH+5, comm));
|
2012-08-20 20:10:23 +06:00
|
|
|
|
|
|
|
PrintMessage ( 3, "now wait ...");
|
2011-07-15 03:36:19 +06:00
|
|
|
|
|
|
|
MPI_Waitall (sendrequests.Size(), &sendrequests[0], MPI_STATUS_IGNORE);
|
2011-07-17 05:13:26 +06:00
|
|
|
|
2019-02-26 15:18:42 +05:00
|
|
|
PrintMessage ( 3, "Sending names");
|
|
|
|
|
|
|
|
sendrequests.SetSize(3*ntasks);
|
|
|
|
/** Send bc/mat/cd*-names **/
|
|
|
|
// nr of names
|
|
|
|
int nnames[4] = {0,0,0,0};
|
|
|
|
nnames[0] = materials.Size();
|
|
|
|
nnames[1] = bcnames.Size();
|
|
|
|
nnames[2] = GetNCD2Names();
|
|
|
|
nnames[3] = GetNCD3Names();
|
|
|
|
int tot_nn = nnames[0] + nnames[1] + nnames[2] + nnames[3];
|
|
|
|
for( int k = 1; k < ntasks; k++)
|
|
|
|
(void) MPI_Isend(nnames, 4, MPI_INT, k, MPI_TAG_MESH+6, comm, &sendrequests[k]);
|
|
|
|
auto iterate_names = [&](auto func) {
|
|
|
|
for (int k = 0; k < nnames[0]; k++) func(materials[k]);
|
|
|
|
for (int k = 0; k < nnames[1]; k++) func(bcnames[k]);
|
|
|
|
for (int k = 0; k < nnames[2]; k++) func(cd2names[k]);
|
|
|
|
for (int k = 0; k < nnames[3]; k++) func(cd3names[k]);
|
|
|
|
};
|
|
|
|
// sizes of names
|
|
|
|
Array<int> name_sizes(tot_nn);
|
|
|
|
tot_nn = 0;
|
|
|
|
iterate_names([&](auto ptr) { name_sizes[tot_nn++] = (ptr==NULL) ? 0 : ptr->size(); });
|
|
|
|
for( int k = 1; k < ntasks; k++)
|
|
|
|
(void) MPI_Isend(&name_sizes[0], tot_nn, MPI_INT, k, MPI_TAG_MESH+6, comm, &sendrequests[ntasks+k]);
|
|
|
|
// names
|
|
|
|
int strs = 0;
|
|
|
|
iterate_names([&](auto ptr) { strs += (ptr==NULL) ? 0 : ptr->size(); });
|
2019-02-28 16:06:30 +05:00
|
|
|
Array<char> compiled_names(strs);
|
2019-02-26 15:18:42 +05:00
|
|
|
strs = 0;
|
|
|
|
iterate_names([&](auto ptr) {
|
|
|
|
if (ptr==NULL) return;
|
|
|
|
auto& name = *ptr;
|
|
|
|
for (int j=0; j < name.size(); j++) compiled_names[strs++] = name[j];
|
|
|
|
});
|
|
|
|
for( int k = 1; k < ntasks; k++)
|
2019-02-28 16:06:30 +05:00
|
|
|
(void) MPI_Isend(&(compiled_names[0]), strs, MPI_CHAR, k, MPI_TAG_MESH+6, comm, &sendrequests[2*ntasks+k]);
|
2017-06-08 20:50:37 +05:00
|
|
|
|
2019-02-26 15:18:42 +05:00
|
|
|
PrintMessage ( 3, "wait for names");
|
2017-06-08 20:50:37 +05:00
|
|
|
|
|
|
|
MPI_Waitall (sendrequests.Size(), &sendrequests[0], MPI_STATUS_IGNORE);
|
|
|
|
|
2019-01-29 20:10:02 +05:00
|
|
|
MPI_Barrier(comm);
|
2019-02-26 15:18:42 +05:00
|
|
|
|
2019-04-15 20:49:45 +05:00
|
|
|
PrintMessage( 3, "Clean up local memory");
|
|
|
|
|
|
|
|
auto & self = const_cast<Mesh&>(*this);
|
|
|
|
self.points = T_POINTS(0);
|
|
|
|
self.surfelements = T_SURFELEMENTS(0);
|
|
|
|
self.volelements = T_VOLELEMENTS(0);
|
|
|
|
self.segments = Array<Segment, 0, size_t>(0);
|
|
|
|
self.lockedpoints = Array<PointIndex>(0);
|
|
|
|
auto cleanup_ptr = [](auto & ptr) {
|
|
|
|
if (ptr != nullptr) {
|
|
|
|
delete ptr;
|
|
|
|
ptr = nullptr;
|
|
|
|
}
|
|
|
|
};
|
|
|
|
cleanup_ptr(self.boundaryedges);
|
|
|
|
cleanup_ptr(self.boundaryedges);
|
|
|
|
cleanup_ptr(self.segmentht);
|
|
|
|
cleanup_ptr(self.surfelementht);
|
|
|
|
self.openelements = Array<Element2d>(0);
|
|
|
|
self.opensegments = Array<Segment>(0);
|
|
|
|
self.numvertices = 0;
|
|
|
|
self.mlbetweennodes = Array<PointIndices<2>,PointIndex::BASE> (0);
|
|
|
|
self.mlparentelement = Array<int>(0);
|
|
|
|
self.mlparentsurfaceelement = Array<int>(0);
|
|
|
|
self.topology.Update();
|
|
|
|
self.curvedelems = new CurvedElements (self);
|
|
|
|
self.clusters = new AnisotropicClusters (self);
|
|
|
|
self.ident = new Identifications (self);
|
|
|
|
self.BuildElementSearchTree();
|
|
|
|
|
|
|
|
// const_cast<Mesh&>(*this).DeleteMesh();
|
|
|
|
|
2019-02-26 15:18:42 +05:00
|
|
|
PrintMessage( 3, "send mesh complete");
|
2011-07-15 03:36:19 +06:00
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
// slaves receive the mesh from the master
|
|
|
|
void Mesh :: ReceiveParallelMesh ( )
|
|
|
|
{
|
|
|
|
int timer = NgProfiler::CreateTimer ("ReceiveParallelMesh");
|
|
|
|
int timer_pts = NgProfiler::CreateTimer ("Receive points");
|
|
|
|
int timer_els = NgProfiler::CreateTimer ("Receive elements");
|
|
|
|
int timer_sels = NgProfiler::CreateTimer ("Receive surface elements");
|
|
|
|
NgProfiler::RegionTimer reg(timer);
|
|
|
|
|
2019-02-13 02:11:55 +05:00
|
|
|
NgMPI_Comm comm = GetCommunicator();
|
|
|
|
int id = comm.Rank();
|
|
|
|
int ntasks = comm.Size();
|
2019-01-29 20:10:02 +05:00
|
|
|
|
2018-07-19 20:33:26 +05:00
|
|
|
int dim;
|
2019-01-29 20:10:02 +05:00
|
|
|
MyMPI_Bcast(dim, comm);
|
2018-07-19 20:33:26 +05:00
|
|
|
SetDimension(dim);
|
|
|
|
|
|
|
|
// Receive number of local elements
|
|
|
|
int nelloc;
|
|
|
|
MPI_Scatter (NULL, 0, MPI_INT,
|
2019-01-29 20:10:02 +05:00
|
|
|
&nelloc, 1, MPI_INT, 0, comm);
|
2018-07-19 20:33:26 +05:00
|
|
|
paralleltop -> SetNE (nelloc);
|
|
|
|
|
2011-07-15 03:36:19 +06:00
|
|
|
// string st;
|
|
|
|
|
|
|
|
// receive vertices
|
|
|
|
NgProfiler::StartTimer (timer_pts);
|
2012-06-13 15:07:11 +06:00
|
|
|
|
2011-07-15 03:36:19 +06:00
|
|
|
Array<int> verts;
|
2019-01-29 20:10:02 +05:00
|
|
|
MyMPI_Recv (verts, 0, MPI_TAG_MESH+1, comm);
|
2018-07-19 20:33:26 +05:00
|
|
|
|
2011-07-15 03:36:19 +06:00
|
|
|
int numvert = verts.Size();
|
|
|
|
paralleltop -> SetNV (numvert);
|
|
|
|
|
|
|
|
// INDEX_CLOSED_HASHTABLE<int> glob2loc_vert_ht (3*numvert+1);
|
2011-07-17 05:13:26 +06:00
|
|
|
INDEX_HASHTABLE<int> glob2loc_vert_ht (3*numvert+1);
|
2011-07-15 03:36:19 +06:00
|
|
|
|
|
|
|
for (int vert = 0; vert < numvert; vert++)
|
|
|
|
{
|
|
|
|
int globvert = verts[vert];
|
|
|
|
paralleltop->SetLoc2Glob_Vert ( vert+1, globvert );
|
|
|
|
glob2loc_vert_ht.Set (globvert, vert+1);
|
|
|
|
}
|
|
|
|
|
|
|
|
for (int i = 0; i < numvert; i++)
|
|
|
|
AddPoint (netgen::Point<3> (0,0,0));
|
|
|
|
|
|
|
|
MPI_Datatype mptype = MeshPoint::MyGetMPIType();
|
|
|
|
MPI_Status status;
|
2019-01-29 20:10:02 +05:00
|
|
|
MPI_Recv( &points[1], numvert, mptype, 0, MPI_TAG_MESH+1, comm, &status);
|
2018-07-19 20:33:26 +05:00
|
|
|
|
|
|
|
Array<int> pp_data;
|
2019-01-29 20:10:02 +05:00
|
|
|
MyMPI_Recv(pp_data, 0, MPI_TAG_MESH+1, comm);
|
2018-07-19 20:33:26 +05:00
|
|
|
|
|
|
|
int maxidentnr = pp_data[0];
|
|
|
|
auto & idents = GetIdentifications();
|
|
|
|
for (int idnr = 1; idnr < maxidentnr+1; idnr++)
|
|
|
|
{
|
|
|
|
|
|
|
|
idents.SetType(idnr, (Identifications::ID_TYPE)pp_data[idnr]);
|
|
|
|
}
|
|
|
|
|
|
|
|
int offset = 2*maxidentnr+1;
|
|
|
|
for(int idnr = 1; idnr < maxidentnr+1; idnr++)
|
|
|
|
{
|
|
|
|
int npairs = pp_data[maxidentnr+idnr];
|
|
|
|
FlatArray<int> pairdata(2*npairs, &pp_data[offset]);
|
|
|
|
offset += 2*npairs;
|
|
|
|
for (int k = 0; k<npairs; k++) {
|
|
|
|
PointIndex loc1 = glob2loc_vert_ht.Get(pairdata[2*k]);
|
|
|
|
PointIndex loc2 = glob2loc_vert_ht.Get(pairdata[2*k+1]);
|
|
|
|
idents.Add(loc1, loc2, idnr);
|
|
|
|
}
|
|
|
|
}
|
2011-07-15 03:36:19 +06:00
|
|
|
|
|
|
|
Array<int> dist_pnums;
|
2019-01-29 20:10:02 +05:00
|
|
|
MyMPI_Recv (dist_pnums, 0, MPI_TAG_MESH+1, comm);
|
2011-07-15 03:36:19 +06:00
|
|
|
|
|
|
|
for (int hi = 0; hi < dist_pnums.Size(); hi += 3)
|
|
|
|
paralleltop ->
|
2012-06-21 20:19:25 +06:00
|
|
|
SetDistantPNum (dist_pnums[hi+1], dist_pnums[hi]); // , dist_pnums[hi+2]);
|
2011-07-15 03:36:19 +06:00
|
|
|
|
|
|
|
NgProfiler::StopTimer (timer_pts);
|
|
|
|
*testout << "got " << numvert << " vertices" << endl;
|
|
|
|
|
|
|
|
{
|
|
|
|
Element el;
|
|
|
|
|
|
|
|
Array<int> elarray;
|
2019-01-29 20:10:02 +05:00
|
|
|
MyMPI_Recv (elarray, 0, MPI_TAG_MESH+2, comm);
|
2011-07-15 03:36:19 +06:00
|
|
|
|
|
|
|
NgProfiler::RegionTimer reg(timer_els);
|
|
|
|
|
|
|
|
for (int ind = 0, elnum = 1; ind < elarray.Size(); elnum++)
|
|
|
|
{
|
|
|
|
paralleltop->SetLoc2Glob_VolEl ( elnum, elarray[ind++]);
|
|
|
|
|
|
|
|
el.SetIndex(elarray[ind++]);
|
|
|
|
el.SetNP(elarray[ind++]);
|
2018-07-20 17:00:34 +05:00
|
|
|
|
2011-07-15 03:36:19 +06:00
|
|
|
for ( int j = 0; j < el.GetNP(); j++)
|
|
|
|
el[j] = glob2loc_vert_ht.Get (elarray[ind++]);
|
|
|
|
|
|
|
|
AddVolumeElement (el);
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
{
|
|
|
|
Array<double> fddata;
|
2019-01-29 20:10:02 +05:00
|
|
|
MyMPI_Recv (fddata, 0, MPI_TAG_MESH+3, comm);
|
2011-07-15 03:36:19 +06:00
|
|
|
for (int i = 0; i < fddata.Size(); i += 6)
|
|
|
|
{
|
2012-06-13 15:07:11 +06:00
|
|
|
int faceind = AddFaceDescriptor
|
|
|
|
(FaceDescriptor(int(fddata[i]), int(fddata[i+1]), int(fddata[i+2]), 0));
|
2011-07-15 03:36:19 +06:00
|
|
|
GetFaceDescriptor(faceind).SetBCProperty (int(fddata[i+3]));
|
|
|
|
GetFaceDescriptor(faceind).domin_singular = fddata[i+4];
|
|
|
|
GetFaceDescriptor(faceind).domout_singular = fddata[i+5];
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
{
|
|
|
|
NgProfiler::RegionTimer reg(timer_sels);
|
|
|
|
Array<int> selbuf;
|
|
|
|
|
2019-01-29 20:10:02 +05:00
|
|
|
MyMPI_Recv ( selbuf, 0, MPI_TAG_MESH+4, comm);
|
2011-07-15 03:36:19 +06:00
|
|
|
|
|
|
|
int ii = 0;
|
|
|
|
int sel = 0;
|
|
|
|
|
|
|
|
int nlocsel = selbuf[ii++];
|
|
|
|
paralleltop -> SetNSE ( nlocsel );
|
|
|
|
|
|
|
|
while (ii < selbuf.Size()-1)
|
|
|
|
{
|
|
|
|
int globsel = selbuf[ii++];
|
|
|
|
int faceind = selbuf[ii++];
|
2014-04-28 13:09:39 +06:00
|
|
|
//bool isghost = selbuf[ii++];
|
2011-07-15 03:36:19 +06:00
|
|
|
int nep = selbuf[ii++];
|
|
|
|
Element2d tri(nep);
|
|
|
|
tri.SetIndex(faceind);
|
|
|
|
for(int j = 1; j <= nep; j++)
|
|
|
|
{
|
|
|
|
tri.PNum(j) = glob2loc_vert_ht.Get (selbuf[ii++]);
|
|
|
|
tri.GeomInfoPi(j).trignum = selbuf[ii++];
|
|
|
|
}
|
|
|
|
paralleltop->SetLoc2Glob_SurfEl ( sel+1, globsel );
|
|
|
|
AddSurfaceElement (tri);
|
|
|
|
sel ++;
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
{
|
|
|
|
Array<double> segmbuf;
|
2019-01-29 20:10:02 +05:00
|
|
|
MyMPI_Recv ( segmbuf, 0, MPI_TAG_MESH+5, comm);
|
2011-07-15 03:36:19 +06:00
|
|
|
|
|
|
|
Segment seg;
|
|
|
|
int globsegi;
|
|
|
|
int ii = 0;
|
|
|
|
int segi = 1;
|
|
|
|
int nsegloc = int ( segmbuf.Size() / 14 ) ;
|
|
|
|
paralleltop -> SetNSegm ( nsegloc );
|
2012-08-20 20:10:23 +06:00
|
|
|
|
2011-07-15 03:36:19 +06:00
|
|
|
while ( ii < segmbuf.Size() )
|
|
|
|
{
|
|
|
|
globsegi = int (segmbuf[ii++]);
|
|
|
|
seg.si = int (segmbuf[ii++]);
|
2018-07-20 17:00:34 +05:00
|
|
|
|
2011-07-15 03:36:19 +06:00
|
|
|
seg.pnums[0] = glob2loc_vert_ht.Get (int(segmbuf[ii++]));
|
|
|
|
seg.pnums[1] = glob2loc_vert_ht.Get (int(segmbuf[ii++]));
|
|
|
|
seg.geominfo[0].trignum = int( segmbuf[ii++] );
|
|
|
|
seg.geominfo[1].trignum = int ( segmbuf[ii++]);
|
|
|
|
seg.surfnr1 = int ( segmbuf[ii++]);
|
|
|
|
seg.surfnr2 = int ( segmbuf[ii++]);
|
|
|
|
seg.edgenr = int ( segmbuf[ii++]);
|
|
|
|
seg.epgeominfo[0].dist = segmbuf[ii++];
|
|
|
|
seg.epgeominfo[1].edgenr = int (segmbuf[ii++]);
|
|
|
|
seg.epgeominfo[1].dist = segmbuf[ii++];
|
|
|
|
|
|
|
|
seg.singedge_left = segmbuf[ii++];
|
|
|
|
seg.singedge_right = segmbuf[ii++];
|
|
|
|
|
|
|
|
seg.epgeominfo[0].edgenr = seg.epgeominfo[1].edgenr;
|
|
|
|
|
|
|
|
seg.domin = seg.surfnr1;
|
|
|
|
seg.domout = seg.surfnr2;
|
|
|
|
if ( seg.pnums[0] >0 && seg.pnums[1] > 0 )
|
|
|
|
{
|
|
|
|
paralleltop-> SetLoc2Glob_Segm ( segi, globsegi );
|
|
|
|
|
|
|
|
AddSegment (seg);
|
|
|
|
segi++;
|
|
|
|
}
|
|
|
|
|
|
|
|
}
|
|
|
|
}
|
2017-06-08 20:50:37 +05:00
|
|
|
|
|
|
|
/** Recv bc-names **/
|
2019-02-26 15:18:42 +05:00
|
|
|
int nnames[4] = {0,0,0,0};
|
|
|
|
MPI_Recv(nnames, 4, MPI_INT, 0, MPI_TAG_MESH+6, comm, MPI_STATUS_IGNORE);
|
|
|
|
materials.SetSize(nnames[0]);
|
|
|
|
bcnames.SetSize(nnames[1]);
|
|
|
|
cd2names.SetSize(nnames[2]);
|
|
|
|
cd3names.SetSize(nnames[3]);
|
|
|
|
|
|
|
|
int tot_nn = nnames[0] + nnames[1] + nnames[2] + nnames[3];
|
|
|
|
Array<int> name_sizes(tot_nn);
|
|
|
|
MPI_Recv(&name_sizes[0], tot_nn, MPI_INT, 0, MPI_TAG_MESH+6, comm, MPI_STATUS_IGNORE);
|
|
|
|
int tot_size = 0;
|
|
|
|
for (int k = 0; k < tot_nn; k++) tot_size += name_sizes[k];
|
2011-07-15 03:36:19 +06:00
|
|
|
|
2019-02-28 16:06:30 +05:00
|
|
|
Array<char> compiled_names(tot_size);
|
|
|
|
MPI_Recv(&(compiled_names[0]), tot_size, MPI_CHAR, 0, MPI_TAG_MESH+6, comm, MPI_STATUS_IGNORE);
|
2019-02-26 15:18:42 +05:00
|
|
|
|
|
|
|
tot_nn = tot_size = 0;
|
|
|
|
auto write_names = [&] (auto & array) {
|
|
|
|
for (int k = 0; k < array.Size(); k++) {
|
|
|
|
int s = name_sizes[tot_nn];
|
|
|
|
array[k] = new string(&compiled_names[tot_size], s);
|
|
|
|
tot_nn++;
|
|
|
|
tot_size += s;
|
|
|
|
}
|
|
|
|
};
|
|
|
|
write_names(materials);
|
|
|
|
write_names(bcnames);
|
|
|
|
write_names(cd2names);
|
|
|
|
write_names(cd3names);
|
2017-06-09 23:21:55 +05:00
|
|
|
|
2019-01-29 20:10:02 +05:00
|
|
|
MPI_Barrier(comm);
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2011-07-15 03:36:19 +06:00
|
|
|
int timerloc = NgProfiler::CreateTimer ("Update local mesh");
|
|
|
|
int timerloc2 = NgProfiler::CreateTimer ("CalcSurfacesOfNode");
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2011-07-15 03:36:19 +06:00
|
|
|
NgProfiler::RegionTimer regloc(timerloc);
|
2012-08-30 19:40:17 +06:00
|
|
|
stringstream str;
|
2012-09-03 15:49:43 +06:00
|
|
|
str << "p" << id << ": got " << GetNE() << " elements and "
|
2017-06-09 23:56:19 +05:00
|
|
|
<< GetNSE() << " surface elements";
|
2019-02-06 18:35:32 +05:00
|
|
|
PrintMessage(2, str.str());
|
|
|
|
// cout << str.str() << endl;
|
2012-08-30 19:40:17 +06:00
|
|
|
// PrintMessage (2, "Got ", GetNE(), " elements and ", GetNSE(), " surface elements");
|
2011-11-01 18:54:07 +06:00
|
|
|
// PrintMessage (2, "Got ", GetNSE(), " surface elements");
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2011-07-15 03:36:19 +06:00
|
|
|
NgProfiler::StartTimer (timerloc2);
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2011-07-15 03:36:19 +06:00
|
|
|
CalcSurfacesOfNode ();
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2011-07-15 03:36:19 +06:00
|
|
|
NgProfiler::StopTimer (timerloc2);
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2017-03-07 14:37:22 +05:00
|
|
|
topology.Update();
|
2011-07-15 03:36:19 +06:00
|
|
|
clusters -> Update();
|
2012-06-21 20:19:25 +06:00
|
|
|
|
|
|
|
// paralleltop -> UpdateCoarseGrid();
|
2011-07-15 03:36:19 +06:00
|
|
|
|
|
|
|
SetNextMajorTimeStamp();
|
|
|
|
}
|
|
|
|
|
2018-07-20 17:00:34 +05:00
|
|
|
|
2009-01-26 00:37:14 +05:00
|
|
|
|
|
|
|
|
2011-07-15 03:36:19 +06:00
|
|
|
|
|
|
|
// distribute the mesh to the slave processors
|
|
|
|
// call it only for the master !
|
|
|
|
void Mesh :: Distribute ()
|
2009-01-26 00:37:14 +05:00
|
|
|
{
|
2019-02-13 02:11:55 +05:00
|
|
|
NgMPI_Comm comm = GetCommunicator();
|
|
|
|
int id = comm.Rank();
|
|
|
|
int ntasks = comm.Size();
|
2012-07-05 17:05:33 +06:00
|
|
|
|
2011-07-17 05:13:26 +06:00
|
|
|
if (id != 0 || ntasks == 1 ) return;
|
2011-07-04 18:29:02 +06:00
|
|
|
|
2011-07-15 03:36:19 +06:00
|
|
|
#ifdef METIS
|
|
|
|
ParallelMetis ();
|
|
|
|
#else
|
|
|
|
for (ElementIndex ei = 0; ei < GetNE(); ei++)
|
|
|
|
(*this)[ei].SetPartition(ntasks * ei/GetNE() + 1);
|
2009-01-26 00:37:14 +05:00
|
|
|
#endif
|
|
|
|
|
2012-07-05 17:05:33 +06:00
|
|
|
/*
|
2011-07-15 03:36:19 +06:00
|
|
|
for (ElementIndex ei = 0; ei < GetNE(); ei++)
|
|
|
|
*testout << "el(" << ei << ") is in part " << (*this)[ei].GetPartition() << endl;
|
|
|
|
for (SurfaceElementIndex ei = 0; ei < GetNSE(); ei++)
|
|
|
|
*testout << "sel(" << int(ei) << ") is in part " << (*this)[ei].GetPartition() << endl;
|
2012-07-05 17:05:33 +06:00
|
|
|
*/
|
|
|
|
|
|
|
|
// MyMPI_SendCmd ("mesh");
|
2011-07-15 03:36:19 +06:00
|
|
|
SendRecvMesh ();
|
|
|
|
}
|
|
|
|
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2012-08-20 20:10:23 +06:00
|
|
|
#ifdef METIS5
|
|
|
|
void Mesh :: ParallelMetis ( )
|
|
|
|
{
|
|
|
|
PrintMessage (3, "call metis 5 ...");
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2012-08-20 20:10:23 +06:00
|
|
|
int timer = NgProfiler::CreateTimer ("Mesh::Partition");
|
|
|
|
NgProfiler::RegionTimer reg(timer);
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2012-08-20 20:10:23 +06:00
|
|
|
idx_t ne = GetNE() + GetNSE() + GetNSeg();
|
|
|
|
idx_t nn = GetNP();
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2012-08-20 20:10:23 +06:00
|
|
|
Array<idx_t> eptr, eind;
|
|
|
|
for (int i = 0; i < GetNE(); i++)
|
|
|
|
{
|
|
|
|
eptr.Append (eind.Size());
|
|
|
|
const Element & el = VolumeElement(i+1);
|
|
|
|
for (int j = 0; j < el.GetNP(); j++)
|
|
|
|
eind.Append (el[j]-1);
|
|
|
|
}
|
|
|
|
for (int i = 0; i < GetNSE(); i++)
|
|
|
|
{
|
|
|
|
eptr.Append (eind.Size());
|
|
|
|
const Element2d & el = SurfaceElement(i+1);
|
|
|
|
for (int j = 0; j < el.GetNP(); j++)
|
|
|
|
eind.Append (el[j]-1);
|
|
|
|
}
|
|
|
|
for (int i = 0; i < GetNSeg(); i++)
|
|
|
|
{
|
|
|
|
eptr.Append (eind.Size());
|
|
|
|
const Segment & el = LineSegment(i+1);
|
2014-04-28 13:09:39 +06:00
|
|
|
eind.Append (el[0]-1);
|
|
|
|
eind.Append (el[1]-1);
|
2012-08-20 20:10:23 +06:00
|
|
|
}
|
|
|
|
eptr.Append (eind.Size());
|
|
|
|
Array<idx_t> epart(ne), npart(nn);
|
|
|
|
|
2019-02-13 02:11:55 +05:00
|
|
|
idxtype nparts = GetCommunicator().Size()-1;
|
2012-08-20 20:10:23 +06:00
|
|
|
|
2019-02-15 18:34:27 +05:00
|
|
|
vol_partition.SetSize(GetNE());
|
|
|
|
surf_partition.SetSize(GetNSE());
|
|
|
|
seg_partition.SetSize(GetNSeg());
|
2019-02-09 02:12:00 +05:00
|
|
|
if (nparts == 1)
|
|
|
|
{
|
|
|
|
for (int i = 0; i < GetNE(); i++)
|
2019-02-15 18:34:27 +05:00
|
|
|
// VolumeElement(i+1).SetPartition(1);
|
|
|
|
vol_partition[i]= 1;
|
2019-02-09 02:12:00 +05:00
|
|
|
for (int i = 0; i < GetNSE(); i++)
|
2019-02-15 18:34:27 +05:00
|
|
|
// SurfaceElement(i+1).SetPartition(1);
|
|
|
|
surf_partition[i] = 1;
|
2019-02-09 02:12:00 +05:00
|
|
|
for (int i = 0; i < GetNSeg(); i++)
|
2019-02-15 18:34:27 +05:00
|
|
|
// LineSegment(i+1).SetPartition(1);
|
|
|
|
seg_partition[i] = 1;
|
2019-02-09 02:12:00 +05:00
|
|
|
}
|
2013-09-02 18:51:53 +06:00
|
|
|
|
2019-02-09 02:12:00 +05:00
|
|
|
else
|
|
|
|
|
|
|
|
{
|
2013-09-02 18:51:53 +06:00
|
|
|
|
2019-02-09 02:12:00 +05:00
|
|
|
idxtype edgecut;
|
|
|
|
|
|
|
|
idxtype ncommon = 3;
|
|
|
|
METIS_PartMeshDual (&ne, &nn, &eptr[0], &eind[0], NULL, NULL, &ncommon, &nparts,
|
|
|
|
NULL, NULL,
|
|
|
|
&edgecut, &epart[0], &npart[0]);
|
|
|
|
|
|
|
|
/*
|
|
|
|
METIS_PartMeshNodal (&ne, &nn, &eptr[0], &eind[0], NULL, NULL, &nparts,
|
|
|
|
NULL, NULL,
|
|
|
|
&edgecut, &epart[0], &npart[0]);
|
|
|
|
*/
|
|
|
|
PrintMessage (3, "metis complete");
|
|
|
|
// cout << "done" << endl;
|
|
|
|
|
|
|
|
for (int i = 0; i < GetNE(); i++)
|
2019-02-15 18:34:27 +05:00
|
|
|
// VolumeElement(i+1).SetPartition(epart[i] + 1);
|
|
|
|
vol_partition[i]= epart[i] + 1;
|
2019-02-09 02:12:00 +05:00
|
|
|
for (int i = 0; i < GetNSE(); i++)
|
2019-02-15 18:34:27 +05:00
|
|
|
// SurfaceElement(i+1).SetPartition(epart[i+GetNE()] + 1);
|
|
|
|
surf_partition[i] = epart[i+GetNE()] + 1;
|
2019-02-09 02:12:00 +05:00
|
|
|
for (int i = 0; i < GetNSeg(); i++)
|
2019-02-15 18:34:27 +05:00
|
|
|
// LineSegment(i+1).SetPartition(epart[i+GetNE()+GetNSE()] + 1);
|
|
|
|
seg_partition[i] = epart[i+GetNE()+GetNSE()] + 1;
|
2019-02-09 02:12:00 +05:00
|
|
|
}
|
|
|
|
|
|
|
|
|
2018-08-01 13:35:26 +05:00
|
|
|
// surface elements attached to volume elements
|
|
|
|
Array<bool, PointIndex::BASE> boundarypoints (GetNP());
|
|
|
|
boundarypoints = false;
|
|
|
|
|
|
|
|
if(GetDimension() == 3)
|
|
|
|
for (SurfaceElementIndex sei = 0; sei < GetNSE(); sei++)
|
|
|
|
{
|
|
|
|
const Element2d & el = (*this)[sei];
|
|
|
|
for (int j = 0; j < el.GetNP(); j++)
|
|
|
|
boundarypoints[el[j]] = true;
|
|
|
|
}
|
|
|
|
else
|
|
|
|
for (SegmentIndex segi = 0; segi < GetNSeg(); segi++)
|
|
|
|
{
|
|
|
|
const Segment & seg = (*this)[segi];
|
|
|
|
for (int j = 0; j < 2; j++)
|
|
|
|
boundarypoints[seg[j]] = true;
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
// Build Pnt2Element table, boundary points only
|
|
|
|
Array<int, PointIndex::BASE> cnt(GetNP());
|
|
|
|
cnt = 0;
|
|
|
|
|
|
|
|
auto loop_els_2d = [&](auto f) {
|
|
|
|
for (SurfaceElementIndex sei = 0; sei < GetNSE(); sei++)
|
|
|
|
{
|
|
|
|
const Element2d & el = (*this)[sei];
|
|
|
|
for (int j = 0; j < el.GetNP(); j++) {
|
|
|
|
f(el[j], sei);
|
|
|
|
}
|
|
|
|
}
|
|
|
|
};
|
|
|
|
auto loop_els_3d = [&](auto f) {
|
|
|
|
for (ElementIndex ei = 0; ei < GetNE(); ei++)
|
|
|
|
{
|
|
|
|
const Element & el = (*this)[ei];
|
|
|
|
for (int j = 0; j < el.GetNP(); j++)
|
|
|
|
f(el[j], ei);
|
|
|
|
}
|
|
|
|
};
|
|
|
|
auto loop_els = [&](auto f)
|
|
|
|
{
|
|
|
|
if (GetDimension() == 3 )
|
|
|
|
loop_els_3d(f);
|
|
|
|
else
|
|
|
|
loop_els_2d(f);
|
|
|
|
};
|
|
|
|
|
|
|
|
|
|
|
|
loop_els([&](auto vertex, int index)
|
|
|
|
{
|
|
|
|
if(boundarypoints[vertex])
|
|
|
|
cnt[vertex]++;
|
|
|
|
});
|
|
|
|
TABLE<int, PointIndex::BASE> pnt2el(cnt);
|
|
|
|
loop_els([&](auto vertex, int index)
|
|
|
|
{
|
|
|
|
if(boundarypoints[vertex])
|
|
|
|
pnt2el.Add(vertex, index);
|
|
|
|
});
|
|
|
|
|
2013-09-02 18:51:53 +06:00
|
|
|
|
|
|
|
if (GetDimension() == 3)
|
|
|
|
{
|
|
|
|
for (SurfaceElementIndex sei = 0; sei < GetNSE(); sei++)
|
|
|
|
{
|
|
|
|
Element2d & sel = (*this)[sei];
|
|
|
|
PointIndex pi1 = sel[0];
|
2018-08-01 13:35:26 +05:00
|
|
|
// FlatArray<ElementIndex> els = pnt2el[pi1];
|
|
|
|
FlatArray<int> els = pnt2el[pi1];
|
2013-09-02 18:51:53 +06:00
|
|
|
|
2019-02-15 18:34:27 +05:00
|
|
|
// sel.SetPartition (-1);
|
|
|
|
surf_partition[sei] = -1;
|
2013-09-02 18:51:53 +06:00
|
|
|
|
|
|
|
for (int j = 0; j < els.Size(); j++)
|
|
|
|
{
|
2018-08-01 13:35:26 +05:00
|
|
|
const Element & el = (*this)[ElementIndex(els[j])];
|
2013-09-02 18:51:53 +06:00
|
|
|
|
|
|
|
bool hasall = true;
|
|
|
|
|
|
|
|
for (int k = 0; k < sel.GetNP(); k++)
|
|
|
|
{
|
|
|
|
bool haspi = false;
|
|
|
|
for (int l = 0; l < el.GetNP(); l++)
|
|
|
|
if (sel[k] == el[l])
|
|
|
|
haspi = true;
|
|
|
|
|
|
|
|
if (!haspi) hasall = false;
|
|
|
|
}
|
|
|
|
|
|
|
|
if (hasall)
|
|
|
|
{
|
2019-02-15 18:34:27 +05:00
|
|
|
// sel.SetPartition (el.GetPartition());
|
|
|
|
surf_partition[sei] = vol_partition[ElementIndex(els[j])];
|
2013-09-02 18:51:53 +06:00
|
|
|
break;
|
|
|
|
}
|
|
|
|
}
|
2019-02-15 18:34:27 +05:00
|
|
|
// if (sel.GetPartition() == -1)
|
|
|
|
if (surf_partition[sei] == -1)
|
2013-09-02 18:51:53 +06:00
|
|
|
cerr << "no volume element found" << endl;
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
for (SegmentIndex si = 0; si < GetNSeg(); si++)
|
|
|
|
{
|
|
|
|
Segment & sel = (*this)[si];
|
|
|
|
PointIndex pi1 = sel[0];
|
2018-08-01 13:35:26 +05:00
|
|
|
FlatArray<int> els = pnt2el[pi1];
|
2013-09-02 18:51:53 +06:00
|
|
|
|
2019-02-15 18:34:27 +05:00
|
|
|
// sel.SetPartition (-1);
|
|
|
|
seg_partition[si] = -1;
|
2013-09-02 18:51:53 +06:00
|
|
|
|
|
|
|
for (int j = 0; j < els.Size(); j++)
|
|
|
|
{
|
2018-08-01 13:35:26 +05:00
|
|
|
const Element & el = (*this)[ElementIndex(els[j])];
|
2013-09-02 18:51:53 +06:00
|
|
|
|
|
|
|
bool haspi[9] = { false }; // max surfnp
|
|
|
|
|
|
|
|
for (int k = 0; k < 2; k++)
|
|
|
|
for (int l = 0; l < el.GetNP(); l++)
|
|
|
|
if (sel[k] == el[l])
|
|
|
|
haspi[k] = true;
|
|
|
|
|
|
|
|
bool hasall = true;
|
|
|
|
for (int k = 0; k < sel.GetNP(); k++)
|
|
|
|
if (!haspi[k]) hasall = false;
|
|
|
|
|
|
|
|
if (hasall)
|
|
|
|
{
|
2019-02-15 18:34:27 +05:00
|
|
|
// sel.SetPartition (el.GetPartition());
|
|
|
|
seg_partition[si] = vol_partition[ElementIndex(els[j])];
|
2013-09-02 18:51:53 +06:00
|
|
|
break;
|
|
|
|
}
|
|
|
|
}
|
2019-02-15 18:34:27 +05:00
|
|
|
// if (sel.GetPartition() == -1)
|
|
|
|
if (seg_partition[si] == -1)
|
2013-09-02 18:51:53 +06:00
|
|
|
cerr << "no volume element found" << endl;
|
|
|
|
}
|
2018-08-01 13:35:26 +05:00
|
|
|
}
|
|
|
|
else
|
|
|
|
{
|
|
|
|
for (SegmentIndex segi = 0; segi < GetNSeg(); segi++)
|
|
|
|
{
|
|
|
|
Segment & seg = (*this)[segi];
|
2019-02-15 18:34:27 +05:00
|
|
|
// seg.SetPartition(-1);
|
|
|
|
seg_partition[segi] = -1;
|
2018-08-01 13:35:26 +05:00
|
|
|
PointIndex pi1 = seg[0];
|
|
|
|
|
|
|
|
FlatArray<int> sels = pnt2el[pi1];
|
|
|
|
for (int j = 0; j < sels.Size(); j++)
|
|
|
|
{
|
|
|
|
SurfaceElementIndex sei = sels[j];
|
|
|
|
Element2d & se = (*this)[sei];
|
|
|
|
bool found = false;
|
|
|
|
for (int l = 0; l < se.GetNP(); l++ && !found)
|
|
|
|
found |= (se[l]==seg[1]);
|
|
|
|
if(found) {
|
2019-02-15 18:34:27 +05:00
|
|
|
// seg.SetPartition(se.GetPartition());
|
|
|
|
seg_partition[segi] = surf_partition[sei];
|
2018-08-01 13:35:26 +05:00
|
|
|
break;
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
2019-02-15 18:34:27 +05:00
|
|
|
// if (seg.GetPartition() == -1) {
|
|
|
|
if (seg_partition[segi] == -1) {
|
2018-08-01 13:35:26 +05:00
|
|
|
cout << endl << "segi: " << segi << endl;
|
|
|
|
cout << "points: " << seg[0] << " " << seg[1] << endl;
|
|
|
|
cout << "surfels: " << endl << sels << endl;
|
|
|
|
throw NgException("no surface element found");
|
|
|
|
}
|
|
|
|
}
|
2013-09-02 18:51:53 +06:00
|
|
|
|
|
|
|
}
|
2012-08-20 20:10:23 +06:00
|
|
|
}
|
2013-09-02 18:51:53 +06:00
|
|
|
|
2012-08-20 20:10:23 +06:00
|
|
|
#endif
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
2012-10-17 16:31:43 +06:00
|
|
|
|
|
|
|
//========================== weights =================================================================
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
// distribute the mesh to the slave processors
|
|
|
|
// call it only for the master !
|
|
|
|
void Mesh :: Distribute (Array<int> & volume_weights , Array<int> & surface_weights, Array<int> & segment_weights)
|
|
|
|
{
|
2019-02-13 02:11:55 +05:00
|
|
|
NgMPI_Comm comm = GetCommunicator();
|
|
|
|
int id = comm.Rank();
|
|
|
|
int ntasks = comm.Size();
|
2012-10-17 16:31:43 +06:00
|
|
|
|
|
|
|
if (id != 0 || ntasks == 1 ) return;
|
|
|
|
|
|
|
|
#ifdef METIS
|
|
|
|
ParallelMetis (volume_weights, surface_weights, segment_weights);
|
|
|
|
#else
|
|
|
|
for (ElementIndex ei = 0; ei < GetNE(); ei++)
|
|
|
|
(*this)[ei].SetPartition(ntasks * ei/GetNE() + 1);
|
|
|
|
#endif
|
|
|
|
|
|
|
|
/*
|
|
|
|
for (ElementIndex ei = 0; ei < GetNE(); ei++)
|
|
|
|
*testout << "el(" << ei << ") is in part " << (*this)[ei].GetPartition() << endl;
|
|
|
|
for (SurfaceElementIndex ei = 0; ei < GetNSE(); ei++)
|
|
|
|
*testout << "sel(" << int(ei) << ") is in part " << (*this)[ei].GetPartition() << endl;
|
|
|
|
*/
|
|
|
|
|
|
|
|
// MyMPI_SendCmd ("mesh");
|
|
|
|
SendRecvMesh ();
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
#ifdef METIS5
|
|
|
|
void Mesh :: ParallelMetis (Array<int> & volume_weights , Array<int> & surface_weights, Array<int> & segment_weights)
|
|
|
|
{
|
|
|
|
PrintMessage (3, "call metis 5 with weights ...");
|
|
|
|
|
|
|
|
// cout << "segment_weights " << segment_weights << endl;
|
|
|
|
// cout << "surface_weights " << surface_weights << endl;
|
|
|
|
// cout << "volume_weights " << volume_weights << endl;
|
|
|
|
|
|
|
|
int timer = NgProfiler::CreateTimer ("Mesh::Partition");
|
|
|
|
NgProfiler::RegionTimer reg(timer);
|
|
|
|
|
|
|
|
idx_t ne = GetNE() + GetNSE() + GetNSeg();
|
|
|
|
idx_t nn = GetNP();
|
|
|
|
|
|
|
|
Array<idx_t> eptr, eind , nwgt;
|
|
|
|
for (int i = 0; i < GetNE(); i++)
|
|
|
|
{
|
|
|
|
eptr.Append (eind.Size());
|
|
|
|
|
|
|
|
const Element & el = VolumeElement(i+1);
|
|
|
|
|
|
|
|
int ind = el.GetIndex();
|
|
|
|
if (volume_weights.Size()<ind)
|
|
|
|
nwgt.Append(0);
|
|
|
|
else
|
|
|
|
nwgt.Append (volume_weights[ind -1]);
|
|
|
|
|
|
|
|
for (int j = 0; j < el.GetNP(); j++)
|
|
|
|
eind.Append (el[j]-1);
|
|
|
|
}
|
|
|
|
for (int i = 0; i < GetNSE(); i++)
|
|
|
|
{
|
|
|
|
eptr.Append (eind.Size());
|
|
|
|
const Element2d & el = SurfaceElement(i+1);
|
|
|
|
|
|
|
|
|
|
|
|
int ind = el.GetIndex();
|
|
|
|
ind = GetFaceDescriptor(ind).BCProperty();
|
|
|
|
if (surface_weights.Size()<ind)
|
|
|
|
nwgt.Append(0);
|
|
|
|
else
|
|
|
|
nwgt.Append (surface_weights[ind -1]);
|
|
|
|
|
|
|
|
|
|
|
|
for (int j = 0; j < el.GetNP(); j++)
|
|
|
|
eind.Append (el[j]-1);
|
|
|
|
}
|
|
|
|
for (int i = 0; i < GetNSeg(); i++)
|
|
|
|
{
|
|
|
|
eptr.Append (eind.Size());
|
|
|
|
|
|
|
|
const Segment & el = LineSegment(i+1);
|
|
|
|
|
|
|
|
int ind = el.si;
|
|
|
|
if (segment_weights.Size()<ind)
|
|
|
|
nwgt.Append(0);
|
|
|
|
else
|
|
|
|
nwgt.Append (segment_weights[ind -1]);
|
|
|
|
|
|
|
|
eind.Append (el[0]);
|
|
|
|
eind.Append (el[1]);
|
|
|
|
}
|
|
|
|
|
|
|
|
eptr.Append (eind.Size());
|
|
|
|
Array<idx_t> epart(ne), npart(nn);
|
|
|
|
|
2019-02-13 02:11:55 +05:00
|
|
|
idxtype nparts = GetCommunicator().Size()-1;
|
2019-02-15 18:34:27 +05:00
|
|
|
vol_partition.SetSize(GetNE());
|
|
|
|
surf_partition.SetSize(GetNSE());
|
|
|
|
seg_partition.SetSize(GetNSeg());
|
|
|
|
|
2019-02-09 02:12:00 +05:00
|
|
|
if (nparts == 1)
|
|
|
|
{
|
|
|
|
for (int i = 0; i < GetNE(); i++)
|
2019-02-15 18:34:27 +05:00
|
|
|
// VolumeElement(i+1).SetPartition(1);
|
|
|
|
vol_partition[i] = 1;
|
2019-02-09 02:12:00 +05:00
|
|
|
for (int i = 0; i < GetNSE(); i++)
|
2019-02-15 18:34:27 +05:00
|
|
|
// SurfaceElement(i+1).SetPartition(1);
|
|
|
|
surf_partition[i] = 1;
|
2019-02-09 02:12:00 +05:00
|
|
|
for (int i = 0; i < GetNSeg(); i++)
|
2019-02-15 18:34:27 +05:00
|
|
|
// LineSegment(i+1).SetPartition(1);
|
|
|
|
seg_partition[i] = 1;
|
2019-02-09 02:12:00 +05:00
|
|
|
return;
|
|
|
|
}
|
|
|
|
|
|
|
|
|
2016-02-16 23:04:08 +05:00
|
|
|
idxtype edgecut;
|
2012-10-17 16:31:43 +06:00
|
|
|
|
|
|
|
|
2016-02-16 23:04:08 +05:00
|
|
|
idxtype ncommon = 3;
|
2012-10-17 16:31:43 +06:00
|
|
|
METIS_PartMeshDual (&ne, &nn, &eptr[0], &eind[0], &nwgt[0], NULL, &ncommon, &nparts,
|
|
|
|
NULL, NULL,
|
|
|
|
&edgecut, &epart[0], &npart[0]);
|
|
|
|
/*
|
|
|
|
METIS_PartMeshNodal (&ne, &nn, &eptr[0], &eind[0], NULL, NULL, &nparts,
|
|
|
|
NULL, NULL,
|
|
|
|
&edgecut, &epart[0], &npart[0]);
|
|
|
|
*/
|
|
|
|
PrintMessage (3, "metis complete");
|
|
|
|
// cout << "done" << endl;
|
|
|
|
|
|
|
|
for (int i = 0; i < GetNE(); i++)
|
2019-02-15 18:34:27 +05:00
|
|
|
// VolumeElement(i+1).SetPartition(epart[i] + 1);
|
|
|
|
vol_partition[i] = epart[i] + 1;
|
2012-10-17 16:31:43 +06:00
|
|
|
for (int i = 0; i < GetNSE(); i++)
|
2019-02-15 18:34:27 +05:00
|
|
|
// SurfaceElement(i+1).SetPartition(epart[i+GetNE()] + 1);
|
|
|
|
surf_partition[i] = epart[i+GetNE()] + 1;
|
2012-10-17 16:31:43 +06:00
|
|
|
for (int i = 0; i < GetNSeg(); i++)
|
2019-02-15 18:34:27 +05:00
|
|
|
// LineSegment(i+1).SetPartition(epart[i+GetNE()+GetNSE()] + 1);
|
|
|
|
seg_partition[i] = epart[i+GetNE()+GetNSE()] + 1;
|
2012-10-17 16:31:43 +06:00
|
|
|
}
|
|
|
|
#endif
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
//===========================================================================================
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
2012-08-20 20:10:23 +06:00
|
|
|
#ifdef METIS4
|
2011-07-15 03:36:19 +06:00
|
|
|
void Mesh :: ParallelMetis ( )
|
|
|
|
{
|
|
|
|
int timer = NgProfiler::CreateTimer ("Mesh::Partition");
|
|
|
|
NgProfiler::RegionTimer reg(timer);
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2011-07-15 03:36:19 +06:00
|
|
|
PrintMessage (3, "Metis called");
|
|
|
|
|
|
|
|
if (GetDimension() == 2)
|
|
|
|
{
|
|
|
|
PartDualHybridMesh2D ( ); // neloc );
|
|
|
|
return;
|
|
|
|
}
|
2009-01-26 00:37:14 +05:00
|
|
|
|
|
|
|
|
2011-11-01 18:54:07 +06:00
|
|
|
idx_t ne = GetNE();
|
|
|
|
idx_t nn = GetNP();
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2011-07-17 05:13:26 +06:00
|
|
|
if (ntasks <= 2 || ne <= 1)
|
2011-07-15 03:36:19 +06:00
|
|
|
{
|
|
|
|
if (ntasks == 1) return;
|
|
|
|
|
|
|
|
for (int i=1; i<=ne; i++)
|
|
|
|
VolumeElement(i).SetPartition(1);
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2011-07-17 05:13:26 +06:00
|
|
|
for (int i=1; i<=GetNSE(); i++)
|
|
|
|
SurfaceElement(i).SetPartition(1);
|
|
|
|
|
2011-07-15 03:36:19 +06:00
|
|
|
return;
|
|
|
|
}
|
2009-01-26 00:37:14 +05:00
|
|
|
|
|
|
|
|
2011-07-15 03:36:19 +06:00
|
|
|
bool uniform_els = true;
|
2011-07-17 05:13:26 +06:00
|
|
|
|
|
|
|
ELEMENT_TYPE elementtype = TET;
|
2012-06-21 20:19:25 +06:00
|
|
|
for (int el = 1; el <= GetNE(); el++)
|
|
|
|
if (VolumeElement(el).GetType() != elementtype)
|
2011-07-15 03:36:19 +06:00
|
|
|
{
|
|
|
|
uniform_els = false;
|
|
|
|
break;
|
|
|
|
}
|
|
|
|
|
2011-07-17 05:13:26 +06:00
|
|
|
|
2011-07-15 03:36:19 +06:00
|
|
|
if (!uniform_els)
|
|
|
|
{
|
2011-08-22 20:51:21 +06:00
|
|
|
PartHybridMesh ();
|
2011-07-15 03:36:19 +06:00
|
|
|
}
|
2011-08-22 20:51:21 +06:00
|
|
|
else
|
|
|
|
{
|
2011-11-01 18:54:07 +06:00
|
|
|
|
|
|
|
// uniform (TET) mesh, JS
|
|
|
|
int npe = VolumeElement(1).GetNP();
|
|
|
|
Array<idxtype> elmnts(ne*npe);
|
|
|
|
|
|
|
|
int etype;
|
|
|
|
if (elementtype == TET)
|
|
|
|
etype = 2;
|
|
|
|
else if (elementtype == HEX)
|
|
|
|
etype = 3;
|
|
|
|
|
2011-07-15 03:36:19 +06:00
|
|
|
|
2012-06-21 20:19:25 +06:00
|
|
|
for (int i=1; i<=ne; i++)
|
|
|
|
for (int j=1; j<=npe; j++)
|
|
|
|
elmnts[(i-1)*npe+(j-1)] = VolumeElement(i).PNum(j)-1;
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2012-06-21 20:19:25 +06:00
|
|
|
int numflag = 0;
|
|
|
|
int nparts = ntasks-1;
|
|
|
|
int ncommon = 3;
|
|
|
|
int edgecut;
|
|
|
|
Array<idxtype> epart(ne), npart(nn);
|
2011-07-15 03:36:19 +06:00
|
|
|
|
2012-06-21 20:19:25 +06:00
|
|
|
// if ( ntasks == 1 )
|
|
|
|
// {
|
|
|
|
// (*this) = *mastermesh;
|
|
|
|
// nparts = 4;
|
|
|
|
// metis :: METIS_PartMeshDual (&ne, &nn, elmnts, &etype, &numflag, &nparts,
|
|
|
|
// &edgecut, epart, npart);
|
|
|
|
// cout << "done" << endl;
|
|
|
|
|
|
|
|
// cout << "edge-cut: " << edgecut << ", balance: " << metis :: ComputeElementBalance(ne, nparts, epart) << endl;
|
|
|
|
|
|
|
|
// for (int i=1; i<=ne; i++)
|
|
|
|
// {
|
|
|
|
// mastermesh->VolumeElement(i).SetPartition(epart[i-1]);
|
|
|
|
// }
|
|
|
|
|
|
|
|
// return;
|
|
|
|
// }
|
|
|
|
|
|
|
|
|
|
|
|
int timermetis = NgProfiler::CreateTimer ("Metis itself");
|
|
|
|
NgProfiler::StartTimer (timermetis);
|
2011-07-15 03:36:19 +06:00
|
|
|
|
2011-11-01 18:54:07 +06:00
|
|
|
#ifdef METIS4
|
2012-06-21 20:19:25 +06:00
|
|
|
cout << "call metis(4)_PartMeshDual ... " << flush;
|
|
|
|
METIS_PartMeshDual (&ne, &nn, &elmnts[0], &etype, &numflag, &nparts,
|
|
|
|
&edgecut, &epart[0], &npart[0]);
|
2011-11-01 18:54:07 +06:00
|
|
|
#else
|
2012-06-21 20:19:25 +06:00
|
|
|
cout << "call metis(5)_PartMeshDual ... " << endl;
|
2012-08-20 20:10:23 +06:00
|
|
|
// idx_t options[METIS_NOPTIONS];
|
2012-06-21 20:19:25 +06:00
|
|
|
|
|
|
|
Array<idx_t> eptr(ne+1);
|
|
|
|
for (int j = 0; j < ne+1; j++)
|
|
|
|
eptr[j] = 4*j;
|
|
|
|
|
|
|
|
METIS_PartMeshDual (&ne, &nn, &eptr[0], &elmnts[0], NULL, NULL, &ncommon, &nparts,
|
|
|
|
NULL, NULL,
|
|
|
|
&edgecut, &epart[0], &npart[0]);
|
2011-11-01 18:54:07 +06:00
|
|
|
#endif
|
2012-06-21 20:19:25 +06:00
|
|
|
|
|
|
|
NgProfiler::StopTimer (timermetis);
|
|
|
|
|
|
|
|
cout << "complete" << endl;
|
2011-11-01 18:54:07 +06:00
|
|
|
#ifdef METIS4
|
2012-06-21 20:19:25 +06:00
|
|
|
cout << "edge-cut: " << edgecut << ", balance: "
|
|
|
|
<< ComputeElementBalance(ne, nparts, &epart[0]) << endl;
|
2011-11-01 18:54:07 +06:00
|
|
|
#endif
|
2012-06-21 20:19:25 +06:00
|
|
|
|
|
|
|
// partition numbering by metis : 0 ... ntasks - 1
|
|
|
|
// we want: 1 ... ntasks
|
|
|
|
for (int i=1; i<=ne; i++)
|
|
|
|
VolumeElement(i).SetPartition(epart[i-1] + 1);
|
2011-08-22 20:51:21 +06:00
|
|
|
}
|
2012-06-21 20:19:25 +06:00
|
|
|
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2011-07-17 05:13:26 +06:00
|
|
|
for (int sei = 1; sei <= GetNSE(); sei++ )
|
|
|
|
{
|
|
|
|
int ei1, ei2;
|
|
|
|
GetTopology().GetSurface2VolumeElement (sei, ei1, ei2);
|
|
|
|
Element2d & sel = SurfaceElement (sei);
|
|
|
|
|
|
|
|
for (int j = 0; j < 2; j++)
|
|
|
|
{
|
|
|
|
int ei = (j == 0) ? ei1 : ei2;
|
|
|
|
if ( ei > 0 && ei <= GetNE() )
|
|
|
|
{
|
|
|
|
sel.SetPartition (VolumeElement(ei).GetPartition());
|
|
|
|
break;
|
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
2011-07-15 03:36:19 +06:00
|
|
|
}
|
|
|
|
#endif
|
2009-01-26 00:37:14 +05:00
|
|
|
|
|
|
|
|
2011-07-17 05:13:26 +06:00
|
|
|
void Mesh :: PartHybridMesh ()
|
2011-07-15 03:36:19 +06:00
|
|
|
{
|
|
|
|
#ifdef METIS
|
|
|
|
int ne = GetNE();
|
|
|
|
|
|
|
|
int nn = GetNP();
|
2017-03-07 14:37:22 +05:00
|
|
|
int nedges = topology.GetNEdges();
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2011-07-15 03:36:19 +06:00
|
|
|
idxtype *xadj, * adjacency, *v_weights = NULL, *e_weights = NULL;
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2011-07-15 03:36:19 +06:00
|
|
|
int weightflag = 0;
|
|
|
|
int numflag = 0;
|
|
|
|
int nparts = ntasks - 1;
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2011-07-15 03:36:19 +06:00
|
|
|
int options[5];
|
|
|
|
options[0] = 0;
|
|
|
|
int edgecut;
|
|
|
|
idxtype * part;
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2011-07-15 03:36:19 +06:00
|
|
|
xadj = new idxtype[nn+1];
|
|
|
|
part = new idxtype[nn];
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2011-07-15 03:36:19 +06:00
|
|
|
Array<int> cnt(nn+1);
|
|
|
|
cnt = 0;
|
|
|
|
|
|
|
|
for ( int edge = 1; edge <= nedges; edge++ )
|
|
|
|
{
|
|
|
|
int v1, v2;
|
2017-03-07 14:37:22 +05:00
|
|
|
topology.GetEdgeVertices ( edge, v1, v2);
|
2011-07-15 03:36:19 +06:00
|
|
|
cnt[v1-1] ++;
|
|
|
|
cnt[v2-1] ++;
|
|
|
|
}
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2011-07-15 03:36:19 +06:00
|
|
|
xadj[0] = 0;
|
|
|
|
for ( int n = 1; n <= nn; n++ )
|
|
|
|
{
|
|
|
|
xadj[n] = idxtype(xadj[n-1] + cnt[n-1]);
|
|
|
|
}
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2011-07-15 03:36:19 +06:00
|
|
|
adjacency = new idxtype[xadj[nn]];
|
|
|
|
cnt = 0;
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2011-07-15 03:36:19 +06:00
|
|
|
for ( int edge = 1; edge <= nedges; edge++ )
|
|
|
|
{
|
|
|
|
int v1, v2;
|
2017-03-07 14:37:22 +05:00
|
|
|
topology.GetEdgeVertices ( edge, v1, v2);
|
2011-07-15 03:36:19 +06:00
|
|
|
adjacency[ xadj[v1-1] + cnt[v1-1] ] = v2-1;
|
|
|
|
adjacency[ xadj[v2-1] + cnt[v2-1] ] = v1-1;
|
|
|
|
cnt[v1-1]++;
|
|
|
|
cnt[v2-1]++;
|
|
|
|
}
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2011-07-15 03:36:19 +06:00
|
|
|
for ( int vert = 0; vert < nn; vert++ )
|
|
|
|
{
|
2016-02-16 23:04:08 +05:00
|
|
|
FlatArray<idxtype> array ( cnt[vert], &adjacency[ xadj[vert] ] );
|
2011-07-15 03:36:19 +06:00
|
|
|
BubbleSort(array);
|
|
|
|
}
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2011-11-01 18:54:07 +06:00
|
|
|
#ifdef METIS4
|
2011-07-15 03:36:19 +06:00
|
|
|
METIS_PartGraphKway ( &nn, xadj, adjacency, v_weights, e_weights, &weightflag,
|
|
|
|
&numflag, &nparts, options, &edgecut, part );
|
2011-11-01 18:54:07 +06:00
|
|
|
#else
|
|
|
|
cout << "currently not supported (metis5), A" << endl;
|
|
|
|
#endif
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2011-07-15 03:36:19 +06:00
|
|
|
Array<int> nodesinpart(ntasks);
|
2019-02-15 18:34:27 +05:00
|
|
|
vol_partition.SetSize(ne);
|
2011-07-15 03:36:19 +06:00
|
|
|
for ( int el = 1; el <= ne; el++ )
|
|
|
|
{
|
|
|
|
Element & volel = VolumeElement(el);
|
|
|
|
nodesinpart = 0;
|
2009-01-26 00:37:14 +05:00
|
|
|
|
|
|
|
|
2011-07-15 03:36:19 +06:00
|
|
|
int el_np = volel.GetNP();
|
|
|
|
int partition = 0;
|
|
|
|
for ( int i = 0; i < el_np; i++ )
|
|
|
|
nodesinpart[ part[volel[i]-1]+1 ] ++;
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2011-07-15 03:36:19 +06:00
|
|
|
for ( int i = 1; i < ntasks; i++ )
|
|
|
|
if ( nodesinpart[i] > nodesinpart[partition] )
|
|
|
|
partition = i;
|
|
|
|
|
2019-02-15 18:34:27 +05:00
|
|
|
// volel.SetPartition(partition);
|
|
|
|
vol_partition[el-1] = partition;
|
2011-07-15 03:36:19 +06:00
|
|
|
}
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2011-07-15 03:36:19 +06:00
|
|
|
delete [] xadj;
|
|
|
|
delete [] part;
|
|
|
|
delete [] adjacency;
|
|
|
|
#else
|
|
|
|
cout << "parthybridmesh not available" << endl;
|
|
|
|
#endif
|
|
|
|
}
|
2009-01-26 00:37:14 +05:00
|
|
|
|
|
|
|
|
2011-07-15 03:36:19 +06:00
|
|
|
void Mesh :: PartDualHybridMesh ( ) // Array<int> & neloc )
|
|
|
|
{
|
|
|
|
#ifdef METIS
|
|
|
|
int ne = GetNE();
|
|
|
|
|
|
|
|
// int nn = GetNP();
|
|
|
|
// int nedges = topology->GetNEdges();
|
2017-03-07 14:37:22 +05:00
|
|
|
int nfaces = topology.GetNFaces();
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2011-07-15 03:36:19 +06:00
|
|
|
idxtype *xadj, * adjacency, *v_weights = NULL, *e_weights = NULL;
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2011-07-15 03:36:19 +06:00
|
|
|
int weightflag = 0;
|
2014-04-28 13:09:39 +06:00
|
|
|
// int numflag = 0;
|
2011-07-15 03:36:19 +06:00
|
|
|
int nparts = ntasks - 1;
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2011-07-15 03:36:19 +06:00
|
|
|
int options[5];
|
|
|
|
options[0] = 0;
|
|
|
|
int edgecut;
|
|
|
|
idxtype * part;
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2011-07-15 03:36:19 +06:00
|
|
|
Array<int, 0> facevolels1(nfaces), facevolels2(nfaces);
|
|
|
|
facevolels1 = -1;
|
|
|
|
facevolels2 = -1;
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2011-07-15 03:36:19 +06:00
|
|
|
Array<int, 0> elfaces;
|
|
|
|
xadj = new idxtype[ne+1];
|
|
|
|
part = new idxtype[ne];
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2011-07-15 03:36:19 +06:00
|
|
|
Array<int, 0> cnt(ne+1);
|
|
|
|
cnt = 0;
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2011-07-15 03:36:19 +06:00
|
|
|
for ( int el=1; el <= ne; el++ )
|
|
|
|
{
|
|
|
|
Element volel = VolumeElement(el);
|
2017-03-07 14:37:22 +05:00
|
|
|
topology.GetElementFaces(el, elfaces);
|
2011-07-15 03:36:19 +06:00
|
|
|
for ( int i = 0; i < elfaces.Size(); i++ )
|
|
|
|
{
|
|
|
|
if ( facevolels1[elfaces[i]-1] == -1 )
|
|
|
|
facevolels1[elfaces[i]-1] = el;
|
|
|
|
else
|
2009-01-26 00:37:14 +05:00
|
|
|
{
|
2011-07-15 03:36:19 +06:00
|
|
|
facevolels2[elfaces[i]-1] = el;
|
|
|
|
cnt[facevolels1[elfaces[i]-1]-1]++;
|
|
|
|
cnt[facevolels2[elfaces[i]-1]-1]++;
|
2009-01-26 00:37:14 +05:00
|
|
|
}
|
2011-07-15 03:36:19 +06:00
|
|
|
}
|
|
|
|
}
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2011-07-15 03:36:19 +06:00
|
|
|
xadj[0] = 0;
|
|
|
|
for ( int n = 1; n <= ne; n++ )
|
|
|
|
{
|
|
|
|
xadj[n] = idxtype(xadj[n-1] + cnt[n-1]);
|
|
|
|
}
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2011-07-15 03:36:19 +06:00
|
|
|
adjacency = new idxtype[xadj[ne]];
|
|
|
|
cnt = 0;
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2011-07-15 03:36:19 +06:00
|
|
|
for ( int face = 1; face <= nfaces; face++ )
|
|
|
|
{
|
|
|
|
int e1, e2;
|
|
|
|
e1 = facevolels1[face-1];
|
|
|
|
e2 = facevolels2[face-1];
|
|
|
|
if ( e2 == -1 ) continue;
|
|
|
|
adjacency[ xadj[e1-1] + cnt[e1-1] ] = e2-1;
|
|
|
|
adjacency[ xadj[e2-1] + cnt[e2-1] ] = e1-1;
|
|
|
|
cnt[e1-1]++;
|
|
|
|
cnt[e2-1]++;
|
|
|
|
}
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2011-07-15 03:36:19 +06:00
|
|
|
for ( int el = 0; el < ne; el++ )
|
|
|
|
{
|
2016-02-16 23:04:08 +05:00
|
|
|
FlatArray<idxtype> array ( cnt[el], &adjacency[ xadj[el] ] );
|
2011-07-15 03:36:19 +06:00
|
|
|
BubbleSort(array);
|
|
|
|
}
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2011-07-15 03:36:19 +06:00
|
|
|
int timermetis = NgProfiler::CreateTimer ("Metis itself");
|
|
|
|
NgProfiler::StartTimer (timermetis);
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2011-11-01 18:54:07 +06:00
|
|
|
#ifdef METIS4
|
2011-07-15 03:36:19 +06:00
|
|
|
METIS_PartGraphKway ( &ne, xadj, adjacency, v_weights, e_weights, &weightflag,
|
|
|
|
&numflag, &nparts, options, &edgecut, part );
|
2011-11-01 18:54:07 +06:00
|
|
|
#else
|
|
|
|
cout << "currently not supported (metis5), B" << endl;
|
|
|
|
#endif
|
|
|
|
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2011-07-15 03:36:19 +06:00
|
|
|
NgProfiler::StopTimer (timermetis);
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2011-07-15 03:36:19 +06:00
|
|
|
Array<int> nodesinpart(ntasks);
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2019-02-15 18:34:27 +05:00
|
|
|
vol_partition.SetSize(ne);
|
2011-07-15 03:36:19 +06:00
|
|
|
for ( int el = 1; el <= ne; el++ )
|
|
|
|
{
|
|
|
|
// Element & volel = VolumeElement(el);
|
|
|
|
nodesinpart = 0;
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2019-02-15 18:34:27 +05:00
|
|
|
// VolumeElement(el).SetPartition(part[el-1 ] + 1);
|
|
|
|
vol_partition[el-1] = part[el-1 ] + 1;
|
2011-07-15 03:36:19 +06:00
|
|
|
}
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2011-07-15 03:36:19 +06:00
|
|
|
/*
|
|
|
|
for ( int i=1; i<=ne; i++)
|
|
|
|
{
|
|
|
|
neloc[ VolumeElement(i).GetPartition() ] ++;
|
|
|
|
}
|
|
|
|
*/
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2011-07-15 03:36:19 +06:00
|
|
|
delete [] xadj;
|
|
|
|
delete [] part;
|
|
|
|
delete [] adjacency;
|
|
|
|
#else
|
|
|
|
cout << "partdualmesh not available" << endl;
|
|
|
|
#endif
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2011-07-15 03:36:19 +06:00
|
|
|
}
|
2009-01-26 00:37:14 +05:00
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
2011-07-15 03:36:19 +06:00
|
|
|
void Mesh :: PartDualHybridMesh2D ( )
|
|
|
|
{
|
|
|
|
#ifdef METIS
|
2016-02-16 23:04:08 +05:00
|
|
|
idxtype ne = GetNSE();
|
2011-07-15 03:36:19 +06:00
|
|
|
int nv = GetNV();
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2011-07-15 03:36:19 +06:00
|
|
|
Array<idxtype> xadj(ne+1);
|
|
|
|
Array<idxtype> adjacency(ne*4);
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2011-07-15 03:36:19 +06:00
|
|
|
// first, build the vertex 2 element table:
|
|
|
|
Array<int, PointIndex::BASE> cnt(nv);
|
|
|
|
cnt = 0;
|
|
|
|
for (SurfaceElementIndex sei = 0; sei < GetNSE(); sei++)
|
|
|
|
for (int j = 0; j < (*this)[sei].GetNP(); j++)
|
|
|
|
cnt[ (*this)[sei][j] ] ++;
|
|
|
|
|
|
|
|
TABLE<SurfaceElementIndex, PointIndex::BASE> vert2els(cnt);
|
|
|
|
for (SurfaceElementIndex sei = 0; sei < GetNSE(); sei++)
|
|
|
|
for (int j = 0; j < (*this)[sei].GetNP(); j++)
|
|
|
|
vert2els.Add ((*this)[sei][j], sei);
|
|
|
|
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2011-07-15 03:36:19 +06:00
|
|
|
// find all neighbour elements
|
|
|
|
int cntnb = 0;
|
|
|
|
Array<int> marks(ne); // to visit each neighbour just once
|
|
|
|
marks = -1;
|
|
|
|
for (SurfaceElementIndex sei = 0; sei < ne; sei++)
|
2009-01-26 00:37:14 +05:00
|
|
|
{
|
2011-07-15 03:36:19 +06:00
|
|
|
xadj[sei] = cntnb;
|
|
|
|
for (int j = 0; j < (*this)[sei].GetNP(); j++)
|
2009-01-26 00:37:14 +05:00
|
|
|
{
|
2011-07-15 03:36:19 +06:00
|
|
|
PointIndex vnr = (*this)[sei][j];
|
|
|
|
|
|
|
|
// all elements with at least one common vertex
|
|
|
|
for (int k = 0; k < vert2els[vnr].Size(); k++)
|
2009-01-26 00:37:14 +05:00
|
|
|
{
|
2011-07-15 03:36:19 +06:00
|
|
|
SurfaceElementIndex sei2 = vert2els[vnr][k];
|
|
|
|
if (sei == sei2) continue;
|
|
|
|
if (marks[sei2] == sei) continue;
|
|
|
|
|
|
|
|
// neighbour, if two common vertices
|
|
|
|
int common = 0;
|
|
|
|
for (int m1 = 0; m1 < (*this)[sei].GetNP(); m1++)
|
|
|
|
for (int m2 = 0; m2 < (*this)[sei2].GetNP(); m2++)
|
|
|
|
if ( (*this)[sei][m1] == (*this)[sei2][m2])
|
|
|
|
common++;
|
|
|
|
|
|
|
|
if (common >= 2)
|
|
|
|
{
|
|
|
|
marks[sei2] = sei; // mark as visited
|
|
|
|
adjacency[cntnb++] = sei2;
|
|
|
|
}
|
2009-01-26 00:37:14 +05:00
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
2011-07-15 03:36:19 +06:00
|
|
|
xadj[ne] = cntnb;
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2011-07-15 03:36:19 +06:00
|
|
|
idxtype *v_weights = NULL, *e_weights = NULL;
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2016-02-16 23:04:08 +05:00
|
|
|
idxtype weightflag = 0;
|
2014-04-28 13:09:39 +06:00
|
|
|
// int numflag = 0;
|
2016-02-16 23:04:08 +05:00
|
|
|
idxtype nparts = ntasks - 1;
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2016-02-16 23:04:08 +05:00
|
|
|
idxtype edgecut;
|
2011-07-15 03:36:19 +06:00
|
|
|
Array<idxtype> part(ne);
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2011-07-15 03:36:19 +06:00
|
|
|
for ( int el = 0; el < ne; el++ )
|
|
|
|
BubbleSort (adjacency.Range (xadj[el], xadj[el+1]));
|
2011-11-01 18:54:07 +06:00
|
|
|
|
|
|
|
#ifdef METIS4
|
2014-04-28 13:09:39 +06:00
|
|
|
int options[5];
|
|
|
|
options[0] = 0;
|
2011-07-15 03:36:19 +06:00
|
|
|
METIS_PartGraphKway ( &ne, &xadj[0], &adjacency[0], v_weights, e_weights, &weightflag,
|
|
|
|
&numflag, &nparts, options, &edgecut, &part[0] );
|
2011-11-01 18:54:07 +06:00
|
|
|
#else
|
|
|
|
idx_t ncon = 1;
|
|
|
|
METIS_PartGraphKway ( &ne, &ncon, &xadj[0], &adjacency[0],
|
|
|
|
v_weights, NULL, e_weights,
|
|
|
|
&nparts,
|
|
|
|
NULL, NULL, NULL,
|
|
|
|
&edgecut, &part[0] );
|
|
|
|
#endif
|
|
|
|
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2019-02-15 18:34:27 +05:00
|
|
|
surf_partition.SetSize(ne);
|
2011-07-15 03:36:19 +06:00
|
|
|
for (SurfaceElementIndex sei = 0; sei < ne; sei++)
|
2019-02-15 18:34:27 +05:00
|
|
|
// (*this) [sei].SetPartition (part[sei]+1);
|
|
|
|
surf_partition[sei] = part[sei]+1;
|
2011-07-15 03:36:19 +06:00
|
|
|
#else
|
|
|
|
cout << "partdualmesh not available" << endl;
|
2011-07-04 18:29:02 +06:00
|
|
|
#endif
|
|
|
|
|
2011-07-15 03:36:19 +06:00
|
|
|
}
|
2011-07-04 18:29:02 +06:00
|
|
|
|
2009-01-26 00:37:14 +05:00
|
|
|
|
2011-07-15 03:36:19 +06:00
|
|
|
|
2009-01-26 00:37:14 +05:00
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
#endif
|