mirror of
https://github.com/NGSolve/netgen.git
synced 2024-11-11 16:49:16 +05:00
weights for metis distribution
This commit is contained in:
parent
2aa04e7dd5
commit
c2a51d61a4
@ -257,7 +257,63 @@ void Ng_LoadMesh (const char * filename)
|
||||
if (ntasks > 1)
|
||||
{
|
||||
// MyMPI_SendCmd ("mesh");
|
||||
mesh -> Distribute();
|
||||
// mesh -> Distribute();
|
||||
|
||||
char * weightsfilename = new char [strlen(filename)+1];
|
||||
strcpy (weightsfilename, filename);
|
||||
weightsfilename[strlen (weightsfilename)-3] = 'w';
|
||||
weightsfilename[strlen (weightsfilename)-2] = 'e';
|
||||
weightsfilename[strlen (weightsfilename)-1] = 'i';
|
||||
|
||||
ifstream weightsfile(weightsfilename);
|
||||
delete [] weightsfilename;
|
||||
|
||||
if (!(weightsfile.good()))
|
||||
mesh -> Distribute();
|
||||
else
|
||||
{
|
||||
char str[20];
|
||||
bool endfile = false;
|
||||
int n, dummy;
|
||||
|
||||
Array<int> segment_weights;
|
||||
Array<int> surface_weights;
|
||||
Array<int> volume_weights;
|
||||
|
||||
while (weightsfile.good() && !endfile)
|
||||
{
|
||||
weightsfile >> str;
|
||||
|
||||
if (strcmp (str, "edgeweights") == 0)
|
||||
{
|
||||
weightsfile >> n;
|
||||
segment_weights.SetSize(n);
|
||||
for (int i=0; i<n; i++)
|
||||
weightsfile >> dummy >> segment_weights[i];
|
||||
}
|
||||
|
||||
if (strcmp (str, "surfaceweights") == 0)
|
||||
{
|
||||
weightsfile >> n;
|
||||
surface_weights.SetSize(n);
|
||||
for (int i=0; i<n; i++)
|
||||
weightsfile >> dummy >> surface_weights[i];
|
||||
}
|
||||
|
||||
if (strcmp (str, "volumeweights") == 0)
|
||||
{
|
||||
weightsfile >> n;
|
||||
volume_weights.SetSize(n);
|
||||
for (int i=0; i<n; i++)
|
||||
weightsfile >> dummy >> volume_weights[i];
|
||||
}
|
||||
|
||||
if (strcmp (str, "endfile") == 0)
|
||||
endfile = true;
|
||||
}
|
||||
|
||||
mesh -> Distribute(volume_weights, surface_weights, segment_weights);
|
||||
}
|
||||
}
|
||||
#endif
|
||||
}
|
||||
|
@ -720,6 +720,8 @@ namespace netgen
|
||||
|
||||
/// distributes the master-mesh to local meshes
|
||||
void Distribute ();
|
||||
void Distribute (Array<int> & volume_weights, Array<int> & surface_weights,
|
||||
Array<int> & segment_weights);
|
||||
|
||||
|
||||
/// find connection to parallel meshes
|
||||
@ -730,6 +732,9 @@ namespace netgen
|
||||
|
||||
/// use metis to decompose master mesh
|
||||
void ParallelMetis (); // Array<int> & neloc );
|
||||
void ParallelMetis (Array<int> & volume_weights, Array<int> & surface_weights,
|
||||
Array<int> & segment_weights);
|
||||
|
||||
void PartHybridMesh (); // Array<int> & neloc );
|
||||
void PartDualHybridMesh (); // Array<int> & neloc );
|
||||
void PartDualHybridMesh2D (); // ( Array<int> & neloc );
|
||||
|
@ -822,6 +822,143 @@ namespace netgen
|
||||
|
||||
|
||||
|
||||
|
||||
//========================== 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)
|
||||
{
|
||||
MPI_Comm_size(MPI_COMM_WORLD, &ntasks);
|
||||
MPI_Comm_rank(MPI_COMM_WORLD, &id);
|
||||
|
||||
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);
|
||||
|
||||
int nparts = ntasks-1;
|
||||
int edgecut;
|
||||
|
||||
|
||||
int ncommon = 3;
|
||||
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++)
|
||||
VolumeElement(i+1).SetPartition(epart[i] + 1);
|
||||
for (int i = 0; i < GetNSE(); i++)
|
||||
SurfaceElement(i+1).SetPartition(epart[i+GetNE()] + 1);
|
||||
for (int i = 0; i < GetNSeg(); i++)
|
||||
LineSegment(i+1).SetPartition(epart[i+GetNE()+GetNSE()] + 1);
|
||||
}
|
||||
#endif
|
||||
|
||||
|
||||
|
||||
//===========================================================================================
|
||||
|
||||
|
||||
|
||||
|
||||
|
||||
|
||||
|
||||
|
||||
|
||||
#ifdef METIS4
|
||||
void Mesh :: ParallelMetis ( )
|
||||
{
|
||||
|
Loading…
Reference in New Issue
Block a user