874 lines
28 KiB
C++
874 lines
28 KiB
C++
#include "Mesh.h"
|
|
#include <stdio.h>
|
|
#include <iostream.h>
|
|
#include <iomanip.h>
|
|
#include <fstream.h>
|
|
#include <strstream.h>
|
|
#include <math.h>
|
|
|
|
//============================================================================
|
|
// Mesh constructor. Here we initialize the mesh using the output of "netgen".
|
|
//============================================================================
|
|
|
|
Mesh::Mesh(char *f_name){
|
|
|
|
level = 0;
|
|
|
|
// Read nodes in Z and triangles in TR (and their number), allocate memory
|
|
// for them and for Atribut)
|
|
ReadNetgenOutput(f_name);
|
|
|
|
// Initialize faces F (using TR and Z) and connection array FCon
|
|
InitializeFaces();
|
|
|
|
// Initialize the Dirichlet bdr Dir (and dim_Dir), Atribut (for the nodes)
|
|
// - to do this we use the faces F (their bdr index in FCon).
|
|
InitializeDirichletBdr();
|
|
|
|
// Initialize Structure V and PN.
|
|
InitializeEdgeStructure();
|
|
|
|
// Initialize the edges E, using V and PN
|
|
// InitializeEdges();
|
|
|
|
// Initialize TRFace using (NF, NTR; TR, F, FCon)
|
|
InitializeTRFace();
|
|
|
|
// Initialize refine edge and tetrahedron type using F & TR
|
|
InitializeRefinementEdges();
|
|
|
|
// Print some mesh information:
|
|
int myrank;
|
|
MPI_Comm_rank(MPI_COMM_WORLD, &myrank);
|
|
if (myrank == 0){
|
|
cout << endl
|
|
<<"+============================================+\n"
|
|
<<"| Structure initialization on Level 0 done |\n"
|
|
<<"+============================================+\n"
|
|
<<"| Number of tetrahedrons = " << setw(9) << NTR <<" |\n"
|
|
<<"| Number of faces = " << setw(9) << NF <<" |\n"
|
|
// <<"| Number of edges = " << setw(9) << NE <<" |\n"
|
|
<<"| Number of nodes = "<<setw(9)<<NN[level]<<" |\n"
|
|
<<"| Number of Dirichlet nodes = "<<setw(9)<<dim_Dir[0]<<" |\n"
|
|
<<"+============================================+\n\n";
|
|
}
|
|
}
|
|
|
|
|
|
//============================================================================
|
|
// return the index of Aij in the one dimensional array PN for level l
|
|
//============================================================================
|
|
int Mesh::ind_ij(int l, int i, int j){
|
|
int k, end = V[l][i+1];
|
|
for(k=V[l][i]; k<end; k++)
|
|
if (PN[l][k] == j)
|
|
return k;
|
|
return -1; // for not found
|
|
}
|
|
|
|
|
|
//============================================================================
|
|
// Print the solution in MCGL format. If the default (proc = -1) all the
|
|
// processors send information to processor 0 which send everything through
|
|
// a socket to the visualization program. If proc != -1 only processor proc
|
|
// visualizes it's part of the mesh.
|
|
//============================================================================
|
|
int call_socket(char *hostname, unsigned short portnum);
|
|
void send_through_socket( int s, char *Buf, int size);
|
|
//============================================================================
|
|
|
|
void Mesh::PrintMCGL(char *f_name, double *solution, int proc){
|
|
int i, myrank, numProcs, s, nn, nf, n_bdr_f = 0;
|
|
char *Buf;
|
|
double scale;
|
|
|
|
ostrstream outs;
|
|
|
|
MPI_Status Status;
|
|
MPI_Comm_rank(MPI_COMM_WORLD, &myrank);
|
|
MPI_Comm_size(MPI_COMM_WORLD, &numProcs);
|
|
|
|
if (proc == -1){
|
|
if (myrank == 0){
|
|
double maximum = 0., minimum = 0.;
|
|
for(i=0;i<NN[level];i++){
|
|
if(solution[i] > maximum)
|
|
maximum=solution[i];
|
|
if(solution[i] < minimum)
|
|
minimum=solution[i];
|
|
}
|
|
if (maximum - minimum < 0.0001)
|
|
scale = 1.;
|
|
else
|
|
scale = 0.5*sqrt(volume())/(maximum - minimum);
|
|
// cout << "Scale = " << scale << "\n";
|
|
|
|
s = call_socket("appleton.llnl.gov", PORTNUM);
|
|
if (s == -1)
|
|
cout << "Couldn't establish socket connection with appleton.llnl.gov\n"
|
|
<< "output written in file " << f_name << endl;
|
|
|
|
outs << "(geometry MC_geometry { : MC_mesh })\n"
|
|
<< "(bbox-draw g1 yes)\n(read geometry {"
|
|
<< "\n appearance { +edge }\n define MC_mesh\n\n"
|
|
<< "# This OFF file was generated by MC.\nOFF\n";
|
|
}
|
|
MPI_Bcast(&scale, 1, MPI_DOUBLE, 0, MPI_COMM_WORLD);
|
|
MPI_Reduce(&NN[level], &nn, 1, MPI_INT, MPI_SUM, 0, MPI_COMM_WORLD);
|
|
|
|
for(i=0; i<NF; i++)
|
|
if (F[i].tetr[1]<0)
|
|
n_bdr_f++;
|
|
MPI_Reduce(&n_bdr_f, &nf, 1, MPI_INT, MPI_SUM, 0, MPI_COMM_WORLD);
|
|
|
|
if (myrank == 0)
|
|
outs << nn << " " << nf << " " << 0 << endl;
|
|
|
|
for(i=0;i<NN[level];i++)
|
|
outs << Z[i].coord[0] << " " << Z[i].coord[1]
|
|
<< " " << Z[i].coord[2] << endl;
|
|
|
|
int length1 = outs.pcount(), length2, dim, shift1, shift2, flag = 1;
|
|
char *buf;
|
|
|
|
shift1 = 0; // This is the shift for subdomain 0
|
|
if (myrank != 0)
|
|
MPI_Recv( &shift1, 1, MPI_INT, MPI_ANY_SOURCE,
|
|
MPI_ANY_TAG, MPI_COMM_WORLD, &Status);
|
|
|
|
if (myrank + 1 < numProcs){
|
|
shift2 = NN[level] + shift1;
|
|
MPI_Send(&shift2, 1, MPI_INT, myrank+1, flag, MPI_COMM_WORLD);
|
|
}
|
|
|
|
MPI_Barrier(MPI_COMM_WORLD);
|
|
|
|
if (myrank == 0)
|
|
for(i=1; i<numProcs; i++){
|
|
MPI_Recv(&dim, 1, MPI_INT, i, MPI_ANY_TAG, MPI_COMM_WORLD, &Status);
|
|
buf = new char[dim+2];
|
|
MPI_Recv(buf, dim, MPI_CHAR, i, MPI_ANY_TAG, MPI_COMM_WORLD, &Status);
|
|
buf[dim] = '\0'; // Simbol marking the end
|
|
outs << buf;
|
|
delete [] buf;
|
|
}
|
|
|
|
for(i=0; i<myrank; i++) rand();
|
|
double c1=((double)rand())/RAND_MAX, c2=((double)rand())/RAND_MAX,
|
|
c3=((double)rand())/RAND_MAX, c4=((double)rand())/RAND_MAX;
|
|
|
|
// double s;
|
|
for(i=0; i<NF; i++){
|
|
// s = (solution[F[i].node[0]] + solution[F[i].node[1]] +
|
|
// solution[F[i].node[2]] )/3.0;
|
|
if (F[i].tetr[1]<0)
|
|
outs << 3 << " "
|
|
<< F[i].node[0]+shift1 << " " << F[i].node[1]+shift1 << " "
|
|
<< F[i].node[2]+shift1
|
|
// << " " << 1-(maximum-s)/(maximum - minimum) << " " << 0.
|
|
// << " " << (maximum-s)/(maximum - minimum) << " " << 0. << endl;
|
|
<< " " << c1 << " " << c2 << " " << c3 << " " << c4 << endl;
|
|
}
|
|
|
|
if (myrank!=0){
|
|
Buf = outs.str(); // freeze the dynamic array
|
|
MPI_Send(&length1, 1, MPI_INT, 0, flag, MPI_COMM_WORLD);
|
|
MPI_Send(Buf, length1, MPI_CHAR, 0, flag, MPI_COMM_WORLD);
|
|
length2 = outs.pcount() - length1;
|
|
MPI_Send(&length2, 1, MPI_INT, 0, flag, MPI_COMM_WORLD);
|
|
MPI_Send(&Buf[length1], length2, MPI_CHAR, 0, flag, MPI_COMM_WORLD);
|
|
}
|
|
else{
|
|
for(i=1; i<numProcs; i++){
|
|
MPI_Recv( &dim, 1, MPI_INT, i, MPI_ANY_TAG, MPI_COMM_WORLD, &Status);
|
|
buf = new char[dim+2];
|
|
MPI_Recv(buf, dim, MPI_CHAR, i, MPI_ANY_TAG, MPI_COMM_WORLD, &Status);
|
|
buf[dim] = '\0';
|
|
outs << buf;
|
|
delete [] buf;
|
|
}
|
|
outs << "\n})\n";
|
|
Buf = outs.str(); // freeze the dynamic array
|
|
|
|
if (s!=-1)
|
|
// Send outs.pcount() chars from Buf through socket s
|
|
send_through_socket( s, Buf, outs.pcount() );
|
|
else{
|
|
ofstream out_stream(f_name);
|
|
out_stream << Buf;
|
|
out_stream.close();
|
|
}
|
|
}
|
|
|
|
delete Buf;
|
|
}
|
|
else
|
|
if (myrank == proc)
|
|
PrintMCGL_sequential(f_name, solution);
|
|
|
|
MPI_Barrier(MPI_COMM_WORLD);
|
|
}
|
|
|
|
|
|
|
|
//============================================================================
|
|
// Print the solution in MCGL format. This is the sequential version.
|
|
// s is the socket port. If !=-1 we visualize using MCGL, otherwise print
|
|
// in file given by f_name.
|
|
//============================================================================
|
|
void Mesh::PrintMCGL_sequential(char *f_name, double *solution){
|
|
int i, n_bdr=0, s;
|
|
char *Buf;
|
|
ostrstream outs;
|
|
|
|
s = call_socket("appleton.llnl.gov", PORTNUM);
|
|
|
|
double maximum = 0., minimum = 0.;
|
|
for(i=0;i<NN[level];i++){
|
|
if(solution[i] > maximum)
|
|
maximum=solution[i];
|
|
if(solution[i] < minimum)
|
|
minimum=solution[i];
|
|
}
|
|
double scale;
|
|
if (maximum - minimum < 0.0001)
|
|
scale = 1.;
|
|
else
|
|
scale = 0.5*sqrt(volume())/(maximum - minimum);
|
|
// cout << "Scale = " << scale << endl;
|
|
|
|
outs << "(geometry MC_geometry { : MC_mesh })\n"
|
|
<< "(bbox-draw g1 yes)\n(read geometry {"
|
|
<< "\n appearance { +edge }\n define MC_mesh\n\n"
|
|
<< "# This OFF file was generated by MC.\nOFF\n";
|
|
|
|
for(i=0; i<NF; i++)
|
|
if (F[i].tetr[1]<0)
|
|
n_bdr++;
|
|
|
|
outs << NN[level] << " " << n_bdr << " " << 0 << endl;
|
|
|
|
for(i=0;i<NN[level];i++)
|
|
outs << Z[i].coord[0] << " " << Z[i].coord[1]
|
|
<< " " << Z[i].coord[2] << endl;
|
|
|
|
double c1=((double)rand())/RAND_MAX, c2=((double)rand())/RAND_MAX,
|
|
c3=((double)rand())/RAND_MAX, c4=((double)rand())/RAND_MAX;
|
|
|
|
// double sol;
|
|
for(i=0; i<NF; i++){
|
|
// sol = (solution[F[i].node[0]] + solution[F[i].node[1]] +
|
|
// solution[F[i].node[2]] )/3.0;
|
|
if (F[i].tetr[1]<0)
|
|
outs << 3 << " "
|
|
<< F[i].node[0] << " " << F[i].node[1] << " " << F[i].node[2]
|
|
// << " " << 1-(maximum-sol)/(maximum - minimum) << " " << 0.
|
|
// << " " << (maximum-sol)/(maximum - minimum) << " " << 0. << endl;
|
|
<< " " << c1 << " " << c2 << " " << c3 << " " << c4 << endl;
|
|
}
|
|
|
|
outs << "\n})\n";
|
|
Buf = outs.str(); // freeze the dynamic array
|
|
|
|
if (s!=-1) // Send outs.pcount() chars from Buf through socket s
|
|
send_through_socket( s, Buf, outs.pcount() );
|
|
else{ // write the same in file given by f_name
|
|
ofstream out_stream(f_name);
|
|
out_stream << Buf;
|
|
out_stream.close();
|
|
}
|
|
|
|
delete Buf;
|
|
}
|
|
|
|
|
|
//============================================================================
|
|
// Print the solution in MCGL format.
|
|
//============================================================================
|
|
void Mesh::PrintSol(char *f_name, double *solution){
|
|
FILE *plot;
|
|
int i;
|
|
|
|
plot=fopen(f_name,"w+");
|
|
if (plot==(FILE *)NULL)
|
|
{cout << "file " << f_name << " is not accessible \n"; exit(1);}
|
|
double maximum = 0., minimum = 0.;
|
|
|
|
for(i=0;i<NN[level];i++){
|
|
if(solution[i] > maximum)
|
|
maximum=solution[i];
|
|
if(solution[i] < minimum)
|
|
minimum=solution[i];
|
|
}
|
|
|
|
double scale = 0.5/(maximum - minimum);
|
|
cout << "Scale = " << scale << endl;
|
|
|
|
fprintf(plot,"(geometry MC_geometry { : MC_mesh })\n");
|
|
fprintf(plot,"(bbox-draw g1 yes)\n(read geometry {");
|
|
fprintf(plot,"\n appearance { +edge }\n define MC_mesh\n\n");
|
|
fprintf(plot,"# This OFF file was generated by MC.\nOFF\n");
|
|
fprintf(plot, "%d %d %d\n", NN[level], NF, 0);
|
|
|
|
for(i=0;i<NN[level];i++)
|
|
fprintf(plot, "%f %f %f\n", Z[i].coord[0],Z[i].coord[1],Z[i].coord[2]);
|
|
|
|
double s;
|
|
for(i=0; i<NF; i++){
|
|
s=(solution[F[i].node[0]]+solution[F[i].node[1]]+
|
|
solution[F[i].node[2]])/3.0;
|
|
fprintf(plot, "%d %d %d %d %4.2f %4.2f %4.2f %4.2f\n",
|
|
3, F[i].node[0], F[i].node[1], F[i].node[2],
|
|
1-(maximum-s)/(maximum - minimum), 0.,
|
|
(maximum-s)/(maximum - minimum), 0.);
|
|
}
|
|
|
|
fprintf(plot,"\n})\n");
|
|
fclose(plot);
|
|
}
|
|
|
|
|
|
//============================================================================
|
|
// Vasilevski output
|
|
//============================================================================
|
|
/*
|
|
void Mesh::PrintMeshMaple(char *f_name){
|
|
FILE *plot;
|
|
FILE *read;
|
|
int i, j, k;
|
|
|
|
plot=fopen(f_name,"w+");
|
|
if (plot==(FILE *)NULL)
|
|
{cout << "file " << f_name << " is not accessible \n";exit(1);}
|
|
|
|
read = fopen("AE_element","r");
|
|
if (read==(FILE *)NULL){
|
|
cout << "file AE_element is not accessible \n";
|
|
exit(1);
|
|
}
|
|
|
|
int *group = new int [NTR];
|
|
int g, el;
|
|
for(i=0; i<NTR; i++){
|
|
fscanf(read,"%d%d", &g, &el);
|
|
group[el] = g;
|
|
}
|
|
fclose(read);
|
|
g++;
|
|
|
|
double cc[g][3];
|
|
int number[g];
|
|
|
|
for(i=0;i<g; i++){
|
|
number[i] = 0;
|
|
cc[i][0] = cc[i][1] = cc[i][2] = 0.;
|
|
}
|
|
|
|
// compute the middle
|
|
for(i=0; i<NTR; i++){
|
|
for(j=0; j<4; j++){ // for the 4 nodes
|
|
for(k=0; k<3; k++) // for the 3 coordinates
|
|
cc[ group[i]][k] += Z[TR[i].node[j]].coord[k];
|
|
|
|
number[ group[i]]++;
|
|
}
|
|
}
|
|
|
|
for(i=0; i<g; i++)
|
|
for(k=0; k<3; k++)
|
|
cc[i][k] /= number[i];
|
|
|
|
real *c0, *c1, *c2, scale = 2.5/3.;
|
|
//double *c0, *c1, *c2, scale = 1/3.;
|
|
for(i=0; i<NTR; i++){ // for all the tetrahedrons
|
|
for(j=0; j<4; j++){ // for the faces
|
|
c0 = Z[F[ TR[i].face[j]].node[0]].GetCoord();
|
|
c1 = Z[F[ TR[i].face[j]].node[1]].GetCoord();
|
|
c2 = Z[F[ TR[i].face[j]].node[2]].GetCoord();
|
|
fprintf(plot,"%7.3f %7.3f %7.3f %7.3f %7.3f %7.3f %7.3f %7.3f %7.3f\n",
|
|
cc[group[i]][0]+(c0[0]-cc[group[i]][0])*scale,
|
|
cc[group[i]][1]+(c0[1]-cc[group[i]][1])*scale,
|
|
cc[group[i]][2]+(c0[2]-cc[group[i]][2])*scale,
|
|
cc[group[i]][0]+(c1[0]-cc[group[i]][0])*scale,
|
|
cc[group[i]][1]+(c1[1]-cc[group[i]][1])*scale,
|
|
cc[group[i]][2]+(c1[2]-cc[group[i]][2])*scale,
|
|
cc[group[i]][0]+(c2[0]-cc[group[i]][0])*scale,
|
|
cc[group[i]][1]+(c2[1]-cc[group[i]][1])*scale,
|
|
cc[group[i]][2]+(c2[2]-cc[group[i]][2])*scale);
|
|
}
|
|
}
|
|
|
|
delete [] group;
|
|
fclose(plot);
|
|
}
|
|
*/
|
|
|
|
//============================================================================
|
|
// The following function prints the connectivity gragh which is used as input
|
|
// for gragh partitioning with metis. The format is:
|
|
// #elements #faces 2
|
|
// weight_element_1 list_of _elements_with_which_it's_connected
|
|
// ...
|
|
// The argument gives weights of the elements - comes from local error estima-
|
|
// mation. If (weight == NULL) there are no weights. The weights have to be
|
|
// integers.
|
|
// On the first line the number of tetrahedrons and internal faces have to be
|
|
// printed, i.e. for now, printing NF is wrong - we have to subtract the num-
|
|
// ber of faces on the boundary, the number that we print at the end.
|
|
//============================================================================
|
|
void Mesh::PrintMetis(double *weight){
|
|
FILE *plot;
|
|
int i, j, *tr, neighbour, n_bdr=0;
|
|
|
|
plot=fopen("Output/metis.input", "w+");
|
|
if (plot==(FILE *)NULL)
|
|
{cout << "The output file can not be opened. Exit.\n"; exit(1);}
|
|
|
|
fprintf(plot,"%d %d 10\n", NTR, NF);
|
|
for(i=0; i<NTR; i++){
|
|
if (weight == NULL)
|
|
fprintf(plot, "%3d ", 1);
|
|
else
|
|
fprintf(plot, "%8.5f ", weight[i]);
|
|
|
|
for(j=0; j<4; j++){ // for the 4 faces
|
|
tr = &(F[ TR[i].face[j] ].tetr[0]);
|
|
if (tr[0] == i) neighbour = tr[1];
|
|
else neighbour = tr[0];
|
|
if (neighbour >= 0)
|
|
fprintf(plot,"%d ", neighbour+1);
|
|
else
|
|
n_bdr++;
|
|
}
|
|
fprintf(plot,"\n");
|
|
}
|
|
|
|
cout << "\n Boundary faces = " << n_bdr << endl;
|
|
fclose(plot);
|
|
}
|
|
|
|
|
|
//============================================================================
|
|
// Vasilevski output
|
|
//============================================================================
|
|
void Mesh::PrintMeshMaple(char *f_name){
|
|
FILE *plot;
|
|
int i, j;
|
|
|
|
plot=fopen(f_name,"w+");
|
|
if (plot==(FILE *)NULL)
|
|
{cout << "file " << f_name << " is not accessible \n"; exit(1);}
|
|
|
|
|
|
//double *c0, *c1, *c2, scale = 2.5/3.;
|
|
real *c0, *c1, *c2, scale = 2/3.;
|
|
scale = 1.;
|
|
real cc[3]; // we use this for tetrah's middle point
|
|
for(i=0; i<NTR; i++){ // for all the tetrahedrons
|
|
GetMiddle(i, cc);
|
|
//if ((TR[i].atribut==1)&&(Z[TR[i].node[0]].coord[2]>250.)){
|
|
for(j=0; j<4; j++){ // for the faces
|
|
c0 = Z[F[ TR[i].face[j]].node[0]].GetCoord();
|
|
c1 = Z[F[ TR[i].face[j]].node[1]].GetCoord();
|
|
c2 = Z[F[ TR[i].face[j]].node[2]].GetCoord();
|
|
fprintf(plot,"%7.3f %7.3f %7.3f %7.3f %7.3f %7.3f %7.3f %7.3f %7.3f\n",
|
|
cc[0]+(c0[0]-cc[0])*scale,
|
|
cc[1]+(c0[1]-cc[1])*scale,
|
|
cc[2]+(c0[2]-cc[2])*scale,
|
|
cc[0]+(c1[0]-cc[0])*scale,
|
|
cc[1]+(c1[1]-cc[1])*scale,
|
|
cc[2]+(c1[2]-cc[2])*scale,
|
|
cc[0]+(c2[0]-cc[0])*scale,
|
|
cc[1]+(c2[1]-cc[1])*scale,
|
|
cc[2]+(c2[2]-cc[2])*scale);
|
|
}
|
|
//}
|
|
}
|
|
fclose(plot);
|
|
}
|
|
|
|
|
|
|
|
//============================================================================
|
|
// Print the Mesh.
|
|
//============================================================================
|
|
void Mesh::PrintMesh(char *f_name){
|
|
FILE *plot;
|
|
int i;
|
|
|
|
plot=fopen(f_name,"w+");
|
|
|
|
fprintf(plot,"%d\n", NTR);
|
|
for(i=0; i<NTR; i++)
|
|
fprintf(plot, "%d %5d %5d %5d %5d\n", 1,
|
|
TR[i].node[0]+1, TR[i].node[1]+1, TR[i].node[2]+1,
|
|
TR[i].node[3]+1);
|
|
|
|
fprintf(plot,"%d\n", NN[level]);
|
|
for(i=0;i<NN[level];i++)
|
|
fprintf(plot, "%10.5f %10.5f %10.5f\n",
|
|
Z[i].coord[0],Z[i].coord[1],Z[i].coord[2]);
|
|
|
|
fclose(plot);
|
|
}
|
|
|
|
|
|
//============================================================================
|
|
// Return the length between node i & j
|
|
//============================================================================
|
|
real Mesh::length(int i, int j){
|
|
real res = 0.;
|
|
for(int k=0; k<3; k++)
|
|
res += (Z[i].coord[k]-Z[j].coord[k])*(Z[i].coord[k]-Z[j].coord[k]);
|
|
return sqrt( res);
|
|
}
|
|
|
|
|
|
//============================================================================
|
|
// Print information about tetrahedron i on the screen.
|
|
//============================================================================
|
|
void Mesh::PrintTetrahedron(int k){
|
|
int i, neighbor;
|
|
|
|
cout << "\nTetrahedron " << k << ", type " << TR[k].type << endl;
|
|
printf("Nodes : %4d %4d %4d %4d with coordinates : \n",
|
|
TR[k].node[0],TR[k].node[1], TR[k].node[2], TR[k].node[3]);
|
|
for(i=0; i<4; i++)
|
|
printf("(%4.2f,%4.2f,%4.2f) ", Z[TR[k].node[i]].coord[0],
|
|
Z[TR[k].node[i]].coord[1], Z[TR[k].node[i]].coord[2]);
|
|
printf("\n");
|
|
|
|
for(i=0; i<4; i++){
|
|
neighbor = F[TR[k].face[i]].tetr[0];
|
|
if (neighbor == k) neighbor = F[TR[k].face[i]].tetr[1];
|
|
printf("Faces : %4d - (%4d,%4d,%4d) connection to tetrahedron %3d\n",
|
|
TR[k].face[i], F[TR[k].face[i]].node[0], F[TR[k].face[i]].node[1],
|
|
F[TR[k].face[i]].node[2], neighbor);
|
|
}
|
|
|
|
printf("Refinement face = %d\n", TR[k].face[TR[k].refine]);
|
|
}
|
|
|
|
|
|
//============================================================================
|
|
// Get pointer to the coordinates of node i.
|
|
//============================================================================
|
|
real *Mesh::GetNode(int i){
|
|
return Z[i].coord;
|
|
}
|
|
|
|
//============================================================================
|
|
// Get the node indexes of tetrahedron "i".
|
|
//============================================================================
|
|
int *Mesh::GetTetr(int i){
|
|
return TR[i].node;
|
|
}
|
|
|
|
|
|
//============================================================================
|
|
// Return the volume of tetrahedron i
|
|
//============================================================================
|
|
real Mesh::volume(int i){
|
|
return fabs(TR[i].determinant(Z))/6;
|
|
}
|
|
|
|
//============================================================================
|
|
// Return the coordinates of the center of tetrahedron i.
|
|
//============================================================================
|
|
void Mesh::GetMiddle(int i, real *coord){
|
|
real *c[4];
|
|
int k;
|
|
|
|
for(k=0; k<4; k++)
|
|
c[k] = GetNode(TR[i].node[k]);
|
|
|
|
for(k=0; k<3; k++)
|
|
coord[k] = (c[0][k] + c[1][k] + c[2][k] + c[3][k])*0.25;
|
|
}
|
|
|
|
//============================================================================
|
|
void Mesh::GetMiddleFace(int i, real *coord){
|
|
for(int k=0; k<3; k++)
|
|
coord[k] = ( Z[F[i].node[0]].coord[k] +
|
|
Z[F[i].node[1]].coord[k] +
|
|
Z[F[i].node[2]].coord[k] )/3;
|
|
}
|
|
|
|
//============================================================================
|
|
void Mesh::GetMiddleEdge(int i, int j, real *coord){
|
|
for(int k=0; k<3; k++)
|
|
coord[k] = (Z[i].coord[k] + Z[j].coord[k])/2;
|
|
}
|
|
|
|
//============================================================================
|
|
// Check the conformity of element t.
|
|
//============================================================================
|
|
void Mesh::CheckConformity(int t){
|
|
int i, j, k, found;
|
|
|
|
// Check if the face points belong to tetrahedron t
|
|
for(i=0; i<4; i++) // for the 4 faces
|
|
for(j=0; j<3; j++){ // for the 3 nodes of face i
|
|
found = 0;
|
|
for(k=0; k<4; k++)
|
|
if (F[TR[t].face[i]].node[j]==TR[t].node[k]){
|
|
found = 1;
|
|
break;
|
|
}
|
|
if (!found){
|
|
PrintTetrahedron(t);
|
|
exit(1);
|
|
}
|
|
}
|
|
}
|
|
|
|
|
|
//============================================================================
|
|
// Return the total volume.
|
|
//============================================================================
|
|
real Mesh::volume(){
|
|
int i;
|
|
real vol = 0.;
|
|
|
|
for(i=0; i<NTR; i++)
|
|
vol += volume(i);
|
|
|
|
return vol;
|
|
}
|
|
|
|
//============================================================================
|
|
// Compute the area of face "f". "f" is local (0..3) index of face in tetr.
|
|
//============================================================================
|
|
real Mesh::area(int tetr, int f){
|
|
real *p[4], v[3][3], n[3];
|
|
int i;
|
|
|
|
for(i=1; i<4; i++)
|
|
p[i] = Z[TR[tetr].node[(f+i)%4]].GetCoord();
|
|
|
|
for(i=1; i<3; i++){
|
|
v[i][0] = p[i][0] - p[i+1][0];
|
|
v[i][1] = p[i][1] - p[i+1][1];
|
|
v[i][2] = p[i][2] - p[i+1][2];
|
|
}
|
|
|
|
n[0] = v[1][1]*v[2][2] - v[1][2]*v[2][1];
|
|
n[1] = v[1][2]*v[2][0] - v[1][0]*v[2][2];
|
|
n[2] = v[1][0]*v[2][1] - v[1][1]*v[2][0];
|
|
|
|
return sqrt(n[0]*n[0]+n[1]*n[1]+n[2]*n[2])*0.5;
|
|
}
|
|
|
|
//============================================================================
|
|
// Print the boundary in Maple format.
|
|
//============================================================================
|
|
|
|
void Mesh::PrintBdr(char *f_name){
|
|
FILE *plot;
|
|
int i;
|
|
|
|
plot=fopen(f_name,"w+");
|
|
|
|
if (plot==(FILE *)NULL)
|
|
{cout << "file " << f_name <<" is not accessible \n";exit(1);}
|
|
|
|
for(i=0; i<NF; i++){
|
|
|
|
if (F[i].tetr[1] < 0){
|
|
fprintf(plot,"%7.3f %7.3f %7.3f %7.3f %7.3f %7.3f %7.3f %7.3f %7.3f\n",
|
|
Z[F[i].node[0]].coord[0],Z[F[i].node[0]].coord[1],Z[F[i].node[0]].coord[2],
|
|
Z[F[i].node[1]].coord[0],Z[F[i].node[1]].coord[1],Z[F[i].node[1]].coord[2],
|
|
Z[F[i].node[2]].coord[0],Z[F[i].node[2]].coord[1],Z[F[i].node[2]].coord[2]);
|
|
}
|
|
}
|
|
fclose(plot);
|
|
}
|
|
|
|
//============================================================================
|
|
// Print the solution on the boundary in MCGL format.
|
|
//============================================================================
|
|
void Mesh::PrintSolBdr(char *f_name, double *solution){
|
|
FILE *plot;
|
|
int i;
|
|
|
|
plot=fopen(f_name,"w+");
|
|
if (plot==(FILE *)NULL)
|
|
{cout << "file " << f_name << " is not accessible \n";exit(1);}
|
|
|
|
double maximum = 0., minimum = 0.;
|
|
|
|
for(i=0;i<NN[level];i++){
|
|
if(solution[i] > maximum)
|
|
maximum=solution[i];
|
|
if(solution[i] < minimum)
|
|
minimum=solution[i];
|
|
}
|
|
|
|
double scale = 0.5/(maximum - minimum);
|
|
cout << "Scale = " << scale << endl;
|
|
|
|
fprintf(plot,"(geometry MC_geometry { : MC_mesh })\n");
|
|
fprintf(plot,"(bbox-draw g1 yes)\n(read geometry {");
|
|
fprintf(plot,"\n appearance { +edge }\n define MC_mesh\n\n");
|
|
fprintf(plot,"# This OFF file was generated by MC.\nOFF\n");
|
|
fprintf(plot, "%d %d %d\n", NN[level], NF, 0);
|
|
|
|
for(i=0;i<NN[level];i++)
|
|
fprintf(plot, "%f %f %f\n", Z[i].coord[0],Z[i].coord[1],Z[i].coord[2]);
|
|
|
|
double s;
|
|
int nbdr = 0;
|
|
for(i=0; i<NF; i++){
|
|
if (F[i].tetr[1] < 0){
|
|
nbdr++;
|
|
s=(solution[F[i].node[0]]+solution[F[i].node[1]]+
|
|
solution[F[i].node[2]])/3.0;
|
|
fprintf(plot, "%d %d %d %d %4.2f %4.2f %4.2f %4.2f\n",
|
|
3, F[i].node[0], F[i].node[1], F[i].node[2],
|
|
1-(maximum-s)/(maximum - minimum), 0.,
|
|
(maximum-s)/(maximum - minimum), 0.);
|
|
}
|
|
}
|
|
|
|
fprintf(plot,"\n})\n");
|
|
fclose(plot);
|
|
fprintf(stderr,"Number of boundary faces : %d\n", nbdr);
|
|
}
|
|
|
|
double exact(real *coord);
|
|
//============================================================================
|
|
// This prints a face in Maple format. For example c = 1, d = 10 means to
|
|
// print face y = 10.
|
|
//============================================================================
|
|
void Mesh::PrintFaceMaple(double *solution, char *f_name, int c, double d){
|
|
FILE *plot;
|
|
int i, n1, n2, n3;
|
|
|
|
plot=fopen(f_name,"w+");
|
|
if (plot==(FILE *)NULL)
|
|
{cout << "file " << f_name << " is not accessible \n"; exit(1);}
|
|
|
|
for(i=0; i<NF; i++){
|
|
|
|
if ( (Z[F[i].node[0]].coord[c] == d) &&
|
|
(Z[F[i].node[1]].coord[c] == d) &&
|
|
(Z[F[i].node[2]].coord[c] == d) ){
|
|
/*
|
|
if (F[i].tetr[1] < 0)
|
|
if ( (Z[F[i].node[0]].coord[c] == d) &&
|
|
(Z[F[i].node[1]].coord[c] == d) &&
|
|
(Z[F[i].node[2]].coord[c] == d) ){
|
|
*/
|
|
n1 = F[i].node[0]; n2 = F[i].node[1]; n3 = F[i].node[2];
|
|
fprintf(plot,
|
|
"%7.3f %7.3f %7.3f %7.3f %7.3f %7.3f %7.3f %7.3f %7.3f\n",
|
|
/*
|
|
// This was used to print the error
|
|
Z[n1].coord[(c+2)%3],Z[n1].coord[(c+1)%3],
|
|
fabs(solution[n1]-exact(Z[F[i].node[0]].coord)),
|
|
Z[n2].coord[(c+2)%3],Z[n2].coord[(c+1)%3],
|
|
fabs(solution[n2]-exact(Z[F[i].node[1]].coord)),
|
|
Z[n3].coord[(c+2)%3],Z[n3].coord[(c+1)%3],
|
|
fabs(solution[n3]-exact(Z[F[i].node[2]].coord)));
|
|
*/
|
|
Z[n1].coord[(c+2)%3],Z[n1].coord[(c+1)%3], solution[n1],
|
|
Z[n2].coord[(c+2)%3],Z[n2].coord[(c+1)%3], solution[n2],
|
|
Z[n3].coord[(c+2)%3],Z[n3].coord[(c+1)%3], solution[n3]);
|
|
}
|
|
}
|
|
fclose(plot);
|
|
}
|
|
|
|
|
|
|
|
//============================================================================
|
|
// This function is used by PrintCutMTVPlot. It writes in file given by the
|
|
// pointer plot the intersection of line p1-p2 with the plane z = 0. The
|
|
// interpolated solution at the intersection is also written (MTV format).
|
|
//============================================================================
|
|
void plot_intersection(double p1[3], int i1, double p2[3], int i2, int &ind,
|
|
double sx, double sy, FILE *plot, double *sol){
|
|
double c;
|
|
|
|
if (!((fabs(p1[2])<0.0001)||(fabs(p2[2])<0.0001)))
|
|
if ((p1[2]>0 && p2[2]<0)||(p1[2]<0 && p2[2]>0)){
|
|
c = p1[2]/(p1[2] - p2[2]);
|
|
// the desired projection is achieved by changing x & y coordinates
|
|
// below (+/-) and adding some shift (sx, sy).
|
|
fprintf(plot,"%7.3f %7.3f %9.5f %d\n", p1[0]+c*(p2[0]-p1[0])+sx,
|
|
-p1[1]-c*(p2[1]-p1[1])+sy, sol[i1]+c*(sol[i2]-sol[i1]), ind++);
|
|
}
|
|
}
|
|
|
|
|
|
//============================================================================
|
|
// This prints a Cut through the domain in MTVPlot format. The cutting plane
|
|
// is given by : Ax + By +Cz + D = 0. The solution is interpolated in the
|
|
// intersection of the plane and the domain.
|
|
//============================================================================
|
|
void Mesh::PrintCutMTVPlot(double *solution, char *f_name,
|
|
double A, double B, double C, double D){
|
|
FILE *plot;
|
|
int i, j, k, n[4], ind;
|
|
real *p[4];
|
|
double T[3][3], cp, sp, ct, st, d, q[4][3];
|
|
|
|
double sx = 0., sy = 0.; // shifts in x and y
|
|
|
|
plot=fopen(f_name,"w+");
|
|
if (plot==(FILE *)NULL)
|
|
{cout << "file " << f_name << " is not accessible \n"; exit(1);}
|
|
|
|
fprintf(plot," $ DATA=CONTCURVE Name = MTVPlot\n");
|
|
fprintf(plot," %% contfill\n %% nsteps = 50\n");
|
|
fprintf(plot," %% toplabel = \"3D BIOSCREEN\"\n");
|
|
fprintf(plot," %% xlabel=\"low biodegradation rate\"\n %% ylabel = \"\"\n");
|
|
fprintf(plot," # %% meshplot = true\n");
|
|
|
|
cp = B/sqrt(A*A+B*B);
|
|
sp = A/sqrt(A*A+B*B);
|
|
ct = C/sqrt(A*A+B*B+C*C);
|
|
st = sqrt(A*A+B*B)/sqrt(A*A+B*B+C*C);
|
|
d = D/sqrt(A*A+B*B+C*C);
|
|
|
|
T[0][0] = cp; T[0][1] = -sp; T[0][2] = 0.;
|
|
T[1][0] = ct*sp; T[1][1] = ct*cp; T[1][2] = -st;
|
|
T[2][0] = st*sp; T[2][1] = st*cp; T[2][2] = ct;
|
|
|
|
ind = 1;
|
|
for(i=0; i<NTR; i++){
|
|
for(j=0; j<4; j++){
|
|
n[j] = TR[i].node[j];
|
|
p[j] = Z[n[j]].GetCoord();
|
|
}
|
|
for(j=0; j<4; j++)
|
|
q[j][2] = T[2][0]*p[j][0] + T[2][1]*p[j][1] + T[2][2]*p[j][2] + d;
|
|
|
|
if ((q[0][2]>0 && q[1][2]>0 && q[2][2]>0 && q[3][2]>0) ||
|
|
(q[0][2]<0 && q[1][2]<0 && q[2][2]<0 && q[3][2]<0))
|
|
continue;
|
|
else { // tetrahedron i intersects the cutting plane
|
|
for(j=0; j<4; j++){ // After the transf tetrahedron i is in q
|
|
q[j][0] = T[0][0]*p[j][0] + T[0][1]*p[j][1] + T[0][2]*p[j][2];
|
|
q[j][1] = T[1][0]*p[j][0] + T[1][1]*p[j][1] + T[1][2]*p[j][2];
|
|
}
|
|
|
|
// Now check the intersection of every edge with z = 0
|
|
for(k=0; k<4; k++)
|
|
if (fabs(q[k][2])<0.0001)
|
|
fprintf(plot,"%7.3f %7.3f %9.5f %d\n",
|
|
q[k][0]+sx, q[k][1]+sy, solution[n[k]],ind++);
|
|
plot_intersection(q[0], n[0], q[1], n[1], ind, sx, sy, plot, solution);
|
|
plot_intersection(q[0], n[0], q[2], n[2], ind, sx, sy, plot, solution);
|
|
plot_intersection(q[0], n[0], q[3], n[3], ind, sx, sy, plot, solution);
|
|
plot_intersection(q[1], n[1], q[2], n[2], ind, sx, sy, plot, solution);
|
|
plot_intersection(q[1], n[1], q[3], n[3], ind, sx, sy, plot, solution);
|
|
plot_intersection(q[2], n[2], q[3], n[3], ind, sx, sy, plot, solution);
|
|
fprintf(plot,"\n");
|
|
}
|
|
}
|
|
fclose(plot);
|
|
}
|
|
|
|
//============================================================================
|