36 int main (
int argc,
char *argv[])
39 const char *mesh_file_1 =
"triple-pt-1.mesh";
40 const char *mesh_file_2 =
"triple-pt-2.mesh";
41 const char *sltn_file_1 =
"triple-pt-1.gf";
42 const char *sltn_file_2 =
"triple-pt-2.gf";
43 bool visualization =
true;
48 args.
AddOption(&mesh_file_1,
"-m1",
"--mesh1",
49 "Mesh file for solution 1.");
50 args.
AddOption(&mesh_file_2,
"-m2",
"--mesh2",
51 "Mesh file for solution 2.");
52 args.
AddOption(&sltn_file_1,
"-s1",
"--solution1",
53 "Grid function for solution 1.");
54 args.
AddOption(&sltn_file_2,
"-s2",
"--solution2",
55 "Grid function for solution 2.");
56 args.
AddOption(&pts_cnt_1D,
"-p",
"--points1D",
57 "Number of comparison points in one direction");
58 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
60 "Enable or disable GLVis visualization.");
70 Mesh mesh_1(mesh_file_1, 1, 1,
false);
71 Mesh mesh_2(mesh_file_2, 1, 1,
false);
74 MFEM_VERIFY(dim > 1,
"GSLIB requires a 2D or a 3D mesh" );
77 const int mesh_poly_deg = mesh_1.
GetNodes()->FESpace()->GetOrder(0);
78 cout <<
"Mesh curvature: "
79 << mesh_1.
GetNodes()->OwnFEC()->Name() <<
" " << mesh_poly_deg << endl;
83 MFEM_VERIFY(mesh_poly_deg > 0,
"The order of the mesh must be positive.");
85 cout <<
"Generating equidistant points for:\n"
86 <<
" x in [" << pos_min(0) <<
", " << pos_max(0) <<
"]\n"
87 <<
" y in [" << pos_min(1) <<
", " << pos_max(1) <<
"]\n";
90 cout <<
"z in [" << pos_min(2) <<
", " << pos_max(2) <<
"]\n";
93 ifstream mat_stream_1(sltn_file_1), mat_stream_2(sltn_file_2);
100 char vishost[] =
"localhost";
103 sout1.
open(vishost, visport);
104 sout2.
open(vishost, visport);
107 cout <<
"Unable to connect to GLVis server at "
108 << vishost <<
':' << visport << endl;
113 sout1 <<
"solution\n" << mesh_1 << func_1
114 <<
"window_title 'Solution 1'"
115 <<
"window_geometry 0 0 600 600";
116 if (dim == 2) { sout1 <<
"keys RmjAc"; }
117 if (dim == 3) { sout1 <<
"keys mA\n"; }
121 sout2 <<
"solution\n" << mesh_2 << func_2
122 <<
"window_title 'Solution 2'"
123 <<
"window_geometry 600 0 600 600";
124 if (dim == 2) { sout2 <<
"keys RmjAc"; }
125 if (dim == 3) { sout2 <<
"keys mA\n"; }
133 const int pts_cnt = pow(pts_cnt_1D, dim);
134 Vector vxyz(pts_cnt * dim);
142 vxyz(i) = pos_min(0) + ip.
x * (pos_max(0)-pos_min(0));
143 vxyz(pts_cnt + i) = pos_min(1) + ip.
y * (pos_max(1)-pos_min(1));
153 vxyz(i) = pos_min(0) + ip.
x * (pos_max(0)-pos_min(0));
154 vxyz(pts_cnt + i) = pos_min(1) + ip.
y * (pos_max(1)-pos_min(1));
155 vxyz(2*pts_cnt + i) = pos_min(2) + ip.
z * (pos_max(2)-pos_min(2));
160 const double rel_bbox_el = 0.05;
161 const double newton_tol = 1.0e-12;
162 const int npts_at_once = 256;
164 Vector pos_r_out(pts_cnt * dim), dist_p_out(pts_cnt);
165 Vector interp_vals_1(pts_cnt), interp_vals_2(pts_cnt);
168 finder.
Setup(mesh_1, rel_bbox_el, newton_tol, npts_at_once);
169 finder.
FindPoints(vxyz, code_out, task_id_out,
170 el_id_out, pos_r_out, dist_p_out);
171 finder.
Interpolate(code_out, task_id_out, el_id_out,
172 pos_r_out, func_1, interp_vals_1);
175 finder.
Setup(mesh_2, rel_bbox_el, newton_tol, npts_at_once);
176 finder.
FindPoints(vxyz, code_out, task_id_out,
177 el_id_out, pos_r_out, dist_p_out);
178 finder.
Interpolate(code_out, task_id_out, el_id_out,
179 pos_r_out, func_2, interp_vals_2);
182 double avg_diff = 0.0, max_diff = 0.0, diff_p;
183 for (
int p = 0; p < pts_cnt; p++)
185 diff_p = fabs(interp_vals_1(p) - interp_vals_2(p));
187 if (diff_p > max_diff) { max_diff = diff_p; }
192 double *nd1 = n1->
GetData(), *nd2 = n2->GetData();
193 double avg_dist = 0.0;
194 const int node_cnt = n1->
Size() /
dim;
195 if (n1->
Size() == n2->Size())
197 for (
int i = 0; i < node_cnt; i++)
200 for (
int d = 0; d <
dim; d++)
202 const int j = i + d * node_cnt;
203 diff_i += (nd1[j] - nd2[j]) * (nd1[j] - nd2[j]);
205 avg_dist += sqrt(diff_i);
207 avg_dist /= node_cnt;
209 else { avg_dist = -1.0; }
211 std::cout <<
"Avg position difference: " << avg_dist << std::endl
212 <<
"Searched " << pts_cnt <<
" points.\n"
213 <<
"Max diff: " << max_diff << std::endl
214 <<
"Avg diff: " << avg_diff << std::endl;
220 const int nodes_cnt = vxyz.Size() /
dim;
223 el_id_out.SetSize(nodes_cnt); code_out.SetSize(nodes_cnt);
224 task_id_out.
SetSize(nodes_cnt);
225 pos_r_out.SetSize(nodes_cnt * dim); dist_p_out.
SetSize(nodes_cnt * dim);
226 interp_vals_2.
SetSize(nodes_cnt);
227 finder.
Setup(mesh_2, rel_bbox_el, newton_tol, npts_at_once);
228 finder.
FindPoints(vxyz, code_out, task_id_out,
229 el_id_out, pos_r_out, dist_p_out);
230 finder.
Interpolate(code_out, task_id_out, el_id_out,
231 pos_r_out, func_2, interp_vals_2);
232 for (
int n = 0; n < nodes_cnt; n++)
234 diff(n) = fabs(func_1(n) - interp_vals_2(n));
239 char vishost[] =
"localhost";
242 sout3.
open(vishost, visport);
244 sout3 <<
"solution\n" << mesh_1 << diff
245 <<
"window_title 'Difference'"
246 <<
"window_geometry 1200 0 600 600";
247 if (dim == 2) { sout3 <<
"keys RmjAcpppppppppppppppppppppp"; }
248 if (dim == 3) { sout3 <<
"keys mA\n"; }
257 const double vol_diff = diff * lf;
258 std::cout <<
"Vol diff: " << vol_diff << std::endl;
int GetNPoints() const
Returns the number of the points in the integration rule.
void FindPoints(Vector &point_pos, Array< unsigned int > &codes, Array< unsigned int > &proc_ids, Array< unsigned int > &elem_ids, Vector &ref_pos, Vector &dist)
Class for domain integration L(v) := (f, v)
Class for an integration rule - an Array of IntegrationPoint.
Class for grid function - Vector with associated FE space.
Subclass constant coefficient.
void SetSize(int s)
Resize the vector to size s.
void Setup(Mesh &m, double bb_t, double newt_tol, int npt_max)
void GetBoundingBox(Vector &min, Vector &max, int ref=2)
Returns the minimum and maximum corners of the mesh bounding box.
int Size() const
Returns the size of the vector.
int main(int argc, char *argv[])
double * GetData() const
Return a pointer to the beginning of the Vector data.
void Interpolate(Array< unsigned int > &codes, Array< unsigned int > &proc_ids, Array< unsigned int > &elem_ids, Vector &ref_pos, const GridFunction &field_in, Vector &field_out)
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
virtual void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1)
const IntegrationRule & GetNodes() const
Nodes: x_i = i/(n-1), i=0,...,n-1.
FiniteElementSpace * FESpace()
void PrintUsage(std::ostream &out) const
void AddOption(bool *var, const char *enable_short_name, const char *enable_long_name, const char *disable_short_name, const char *disable_long_name, const char *description, bool required=false)
void SetSize(int nsize)
Change logical size of the array, keep existing entries.
Class for integration point with weight.
void PrintOptions(std::ostream &out) const
int open(const char hostname[], int port)
void GetNodes(Vector &node_coord) const