44 for (
int d = 0; d <
dim; d++) { res += x(d) * x(d); }
51 for (
int i = 1; i < F.
Size(); i++) { F(i) = (i+1)*pow(-1, i)*F(0); }
54 int main (
int argc,
char *argv[])
57 const char *src_mesh_file =
"../meshing/square01.mesh";
58 const char *tar_mesh_file =
"../../data/inline-tri.mesh";
59 const char *src_sltn_file =
"must_be_provided_by_the_user.gf";
60 int src_fieldtype = 0;
65 bool visualization =
true;
69 args.
AddOption(&src_mesh_file,
"-m1",
"--mesh1",
70 "Mesh file for the starting solution.");
71 args.
AddOption(&tar_mesh_file,
"-m2",
"--mesh2",
72 "Mesh file for interpolation.");
73 args.
AddOption(&src_sltn_file,
"-s1",
"--solution1",
74 "(optional) GridFunction file compatible with src_mesh_file."
75 "Set src_fieldtype to -1 if this option is used.");
76 args.
AddOption(&src_fieldtype,
"-fts",
"--field-type-src",
77 "Source GridFunction type:"
78 "0 - H1 (default), 1 - L2, 2 - H(div), 3 - H(curl).");
79 args.
AddOption(&src_ncomp,
"-nc",
"--ncomp",
80 "Number of components for H1 or L2 GridFunctions.");
81 args.
AddOption(&ref_levels,
"-r",
"--refine",
82 "Number of refinements of the interpolation mesh.");
83 args.
AddOption(&fieldtype,
"-ft",
"--field-type",
84 "Target GridFunction type: -1 - source GridFunction type (default),"
85 "0 - H1, 1 - L2, 2 - H(div), 3 - H(curl).");
87 "Order of the interpolated solution.");
88 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
90 "Enable or disable GLVis visualization.");
100 Mesh mesh_1(src_mesh_file, 1, 1,
false);
101 Mesh mesh_2(tar_mesh_file, 1, 1,
false);
103 MFEM_ASSERT(dim == mesh_2.
Dimension(),
"Source and target meshes "
104 "must be in the same dimension.");
105 MFEM_VERIFY(dim > 1,
"GSLIB requires a 2D or a 3D mesh" );
107 for (
int lev = 0; lev < ref_levels; lev++)
114 const int mesh_poly_deg =
115 mesh_2.
GetNodes()->FESpace()->GetElementOrder(0);
116 cout <<
"Source mesh curvature: "
117 << mesh_1.
GetNodes()->OwnFEC()->Name() << endl
118 <<
"Target mesh curvature: "
119 << mesh_2.
GetNodes()->OwnFEC()->Name() << endl;
121 int src_vdim = src_ncomp;
125 if (src_fieldtype < 0)
127 ifstream mat_stream_1(src_sltn_file);
131 else if (src_fieldtype == 0)
135 else if (src_fieldtype == 1)
139 else if (src_fieldtype == 2)
145 else if (src_fieldtype == 3)
153 MFEM_ABORT(
"Invalid FECollection type.");
156 if (src_fieldtype > -1)
171 sout1.
open(vishost, visport);
174 cout <<
"Unable to connect to GLVis server at "
175 << vishost <<
':' << visport << endl;
180 sout1 <<
"solution\n" << mesh_1 << *func_source
181 <<
"window_title 'Source mesh and solution'"
182 <<
"window_geometry 0 0 600 600";
183 if (dim == 2) { sout1 <<
"keys RmjAc"; }
184 if (dim == 3) { sout1 <<
"keys mA\n"; }
193 "are not currently supported.");
197 std::cout <<
"Source FE collection: " << fec_in->
Name() << std::endl;
199 if (src_fieldtype < 0)
205 if (fec_h1) { src_fieldtype = 0; }
206 else if (fec_l2) { src_fieldtype = 1; }
207 else if (fec_rt) { src_fieldtype = 2; }
208 else if (fec_nd) { src_fieldtype = 3; }
209 else { MFEM_ABORT(
"GridFunction type not supported yet."); }
211 if (fieldtype < 0) { fieldtype = src_fieldtype; }
217 int tar_vdim = src_vdim;
221 tar_vdim = (src_fieldtype > 1) ? dim : src_vdim;
223 else if (fieldtype == 1)
226 tar_vdim = (src_fieldtype > 1) ? dim : src_vdim;
228 else if (fieldtype == 2)
232 MFEM_VERIFY(src_fieldtype > 1,
"Cannot interpolate a scalar "
233 "grid function to a vector");
236 else if (fieldtype == 3)
240 MFEM_VERIFY(src_fieldtype > 1,
"Cannot interpolate a scalar "
241 "grid function to a vector");
245 MFEM_ABORT(
"GridFunction type not supported.");
247 std::cout <<
"Target FE collection: " << tar_fec->
Name() << std::endl;
251 const int NE = mesh_2.
GetNE(),
257 if (fieldtype == 0 && order == mesh_poly_deg)
264 for (
int i = 0; i < NE; i++)
273 rowy(vxyz.
GetData() + i*nsp + NE*nsp, nsp),
281 if (dim == 3) { pos.
GetRow(2, rowz); }
284 const int nodes_cnt = vxyz.
Size() /
dim;
287 Vector interp_vals(nodes_cnt*tar_ncomp);
289 finder.
Setup(mesh_1);
290 finder.
Interpolate(vxyz, *func_source, interp_vals);
295 if ((fieldtype == 0 && order == mesh_poly_deg) || fieldtype == 1)
297 func_target = interp_vals;
303 Vector elem_dof_vals(nsp*tar_ncomp);
305 for (
int i = 0; i < mesh_2.
GetNE(); i++)
309 for (
int j = 0; j < nsp; j++)
311 for (
int d = 0; d < tar_ncomp; d++)
314 elem_dof_vals(j+d*nsp) = interp_vals(d*nsp*NE + i*nsp + j);
325 Vector elem_dof_vals(nsp*tar_ncomp);
327 for (
int i = 0; i < mesh_2.
GetNE(); i++)
331 for (
int j = 0; j < nsp; j++)
333 for (
int d = 0; d < tar_ncomp; d++)
336 elem_dof_vals(j*tar_ncomp+d) = interp_vals(d*nsp*NE + i*nsp + j);
352 sout1.
open(vishost, visport);
355 cout <<
"Unable to connect to GLVis server at "
356 << vishost <<
':' << visport << endl;
361 sout1 <<
"solution\n" << mesh_2 << func_target
362 <<
"window_title 'Target mesh and solution'"
363 <<
"window_geometry 600 0 600 600";
364 if (dim == 2) { sout1 <<
"keys RmjAc"; }
365 if (dim == 3) { sout1 <<
"keys mA\n"; }
371 ostringstream rho_name;
372 rho_name <<
"interpolated.gf";
373 ofstream rho_ofs(rho_name.str().c_str());
374 rho_ofs.precision(8);
375 func_target.
Save(rho_ofs);
int GetNPoints() const
Returns the number of the points in the integration rule.
Abstract class for all finite elements.
virtual void Interpolate(const GridFunction &field_in, Vector &field_out)
int Size() const
Return the logical size of the array.
void SetSubVector(const Array< int > &dofs, const double value)
Set the entries listed in dofs to the given value.
Class for an integration rule - an Array of IntegrationPoint.
Class for grid function - Vector with associated FE space.
int GetNumGeometries(int dim) const
Return the number of geometries of the given dimension present in the mesh.
void SetSize(int s)
Resize the vector to size s.
void GetElementVDofs(int i, Array< int > &vdofs) const
Returns indexes of degrees of freedom in array dofs for i'th element.
Data type dense matrix using column-major storage.
int Size() const
Returns the size of the vector.
int GetNE() const
Returns number of elements.
double * GetData() const
Return a pointer to the beginning of the Vector data.
virtual void ProjectFromNodes(Vector &vc, ElementTransformation &Trans, Vector &dofs) const
Given a vector of values at the finite element nodes and a transformation, compute its projection (ap...
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
Geometry::Type GetGeomType() const
Returns the Geometry::Type of the reference element.
double scalar_func(const Vector &x)
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
void GetRow(int r, Vector &row) const
virtual void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1)
void Setup(Mesh &m, const double bb_t=0.1, const double newt_tol=1.0e-12, const int npt_max=256)
const IntegrationRule & GetNodes() const
Get a const reference to the nodes of the element.
int GetVDim() const
Returns vector dimension.
FiniteElementSpace * FESpace()
void PrintUsage(std::ostream &out) const
Print the usage message.
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
A general vector function coefficient.
ElementTransformation * GetElementTransformation(int i) const
Returns ElementTransformation for the i-th element.
double p(const Vector &x, double t)
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
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)
Add a boolean option and set 'var' to receive the value. Enable/disable tags are used to set the bool...
void SetDataAndSize(double *d, int s)
Set the Vector data and size.
virtual const char * Name() const
void vector_func(const Vector &p, Vector &F)
const FiniteElementSpace * GetNodalFESpace() const
virtual const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th element in t...
void PrintOptions(std::ostream &out) const
Print the options.
virtual void ProjectCoefficient(Coefficient &coeff)
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
int open(const char hostname[], int port)
Open the socket stream on 'port' at 'hostname'.
Arbitrary order H(curl)-conforming Nedelec finite elements.
FindPointsGSLIB can robustly evaluate a GridFunction on an arbitrary collection of points...
const FiniteElementCollection * FEColl() const
void GetNodes(Vector &node_coord) const
Arbitrary order H1-conforming (continuous) finite elements.
Arbitrary order "L2-conforming" discontinuous finite elements.
bool Good() const
Return true if the command line options were parsed successfully.