46 for (
int d = 0; d <
dim; d++) { res += x(d) * x(d); }
53 for (
int i = 1; i < F.
Size(); i++) { F(i) = (i+1)*pow(-1, i)*F(0); }
56int main (
int argc,
char *argv[])
59 const char *src_mesh_file =
"../meshing/square01.mesh";
60 const char *tar_mesh_file =
"../../data/inline-tri.mesh";
61 const char *src_sltn_file =
"must_be_provided_by_the_user.gf";
62 int src_fieldtype = 0;
64 int src_gf_ordering = 0;
68 bool visualization =
true;
72 args.
AddOption(&src_mesh_file,
"-m1",
"--mesh1",
73 "Mesh file for the starting solution.");
74 args.
AddOption(&tar_mesh_file,
"-m2",
"--mesh2",
75 "Mesh file for interpolation.");
76 args.
AddOption(&src_sltn_file,
"-s1",
"--solution1",
77 "(optional) GridFunction file compatible with src_mesh_file."
78 "Set src_fieldtype to -1 if this option is used.");
79 args.
AddOption(&src_fieldtype,
"-fts",
"--field-type-src",
80 "Source GridFunction type:"
81 "0 - H1 (default), 1 - L2, 2 - H(div), 3 - H(curl).");
82 args.
AddOption(&src_ncomp,
"-nc",
"--ncomp",
83 "Number of components for H1 or L2 GridFunctions.");
84 args.
AddOption(&src_gf_ordering,
"-gfo",
"--gfo",
85 "Node ordering: 0 (byNodes) or 1 (byVDim)");
86 args.
AddOption(&ref_levels,
"-r",
"--refine",
87 "Number of refinements of the interpolation mesh.");
88 args.
AddOption(&fieldtype,
"-ft",
"--field-type",
89 "Target GridFunction type: -1 - source GridFunction type (default),"
90 "0 - H1, 1 - L2, 2 - H(div), 3 - H(curl).");
92 "Order of the interpolated solution.");
93 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
95 "Enable or disable GLVis visualization.");
105 if (strcmp(src_sltn_file,
"must_be_provided_by_the_user.gf") != 0)
111 Mesh mesh_1(src_mesh_file, 1, 1,
false);
112 Mesh mesh_2(tar_mesh_file, 1, 1,
false);
114 MFEM_ASSERT(
dim == mesh_2.
Dimension(),
"Source and target meshes "
115 "must be in the same dimension.");
116 MFEM_VERIFY(
dim > 1,
"GSLIB requires a 2D or a 3D mesh" );
118 for (
int lev = 0; lev < ref_levels; lev++)
125 const int mesh_poly_deg =
126 mesh_2.
GetNodes()->FESpace()->GetElementOrder(0);
127 cout <<
"Source mesh curvature: "
128 << mesh_1.
GetNodes()->OwnFEC()->Name() << endl
129 <<
"Target mesh curvature: "
130 << mesh_2.
GetNodes()->OwnFEC()->Name() << endl;
132 int src_vdim = src_ncomp;
136 if (src_fieldtype < 0)
138 ifstream mat_stream_1(src_sltn_file);
141 src_fes = func_source->
FESpace();
143 else if (src_fieldtype == 0)
147 else if (src_fieldtype == 1)
151 else if (src_fieldtype == 2)
157 else if (src_fieldtype == 3)
165 MFEM_ABORT(
"Invalid FECollection type.");
168 if (src_fieldtype > -1)
170 MFEM_VERIFY(src_gf_ordering >= 0 && src_gf_ordering <= 1,
171 " Source grid function ordering must be 0 (byNodes)"
189 cout <<
"Unable to connect to GLVis server at "
195 sout1 <<
"solution\n" << mesh_1 << *func_source
196 <<
"window_title 'Source mesh and solution'"
197 <<
"window_geometry 0 0 600 600";
198 if (
dim == 2) { sout1 <<
"keys RmjAc"; }
199 if (
dim == 3) { sout1 <<
"keys mA\n"; }
208 "are not currently supported.");
212 std::cout <<
"Source FE collection: " << fec_in->
Name() << std::endl;
214 if (src_fieldtype < 0)
220 if (fec_h1) { src_fieldtype = 0; }
221 else if (fec_l2) { src_fieldtype = 1; }
222 else if (fec_rt) { src_fieldtype = 2; }
223 else if (fec_nd) { src_fieldtype = 3; }
224 else { MFEM_ABORT(
"GridFunction type not supported yet."); }
226 if (fieldtype < 0) { fieldtype = src_fieldtype; }
232 int tar_vdim = src_vdim;
236 tar_vdim = (src_fieldtype > 1) ?
dim : src_vdim;
238 else if (fieldtype == 1)
241 tar_vdim = (src_fieldtype > 1) ?
dim : src_vdim;
243 else if (fieldtype == 2)
247 MFEM_VERIFY(src_fieldtype > 1,
"Cannot interpolate a scalar "
248 "grid function to a vector");
251 else if (fieldtype == 3)
255 MFEM_VERIFY(src_fieldtype > 1,
"Cannot interpolate a scalar "
256 "grid function to a vector");
260 MFEM_ABORT(
"GridFunction type not supported.");
262 std::cout <<
"Target FE collection: " << tar_fec->
Name() << std::endl;
267 const int NE = mesh_2.
GetNE(),
274 if (fieldtype == 0 && order == mesh_poly_deg)
277 point_ordering = mesh_2.
GetNodes()->FESpace()->GetOrdering();
282 for (
int i = 0; i < NE; i++)
291 rowy(vxyz.
GetData() + i*nsp + NE*nsp, nsp),
303 const int nodes_cnt = vxyz.
Size() /
dim;
306 Vector interp_vals(nodes_cnt*tar_ncomp);
308 finder.
Setup(mesh_1);
309 finder.
Interpolate(vxyz, *func_source, interp_vals, point_ordering);
314 if ((fieldtype == 0 && order == mesh_poly_deg) || fieldtype == 1)
316 func_target = interp_vals;
322 Vector elem_dof_vals(nsp*tar_ncomp);
324 for (
int i = 0; i < mesh_2.
GetNE(); i++)
328 for (
int j = 0; j < nsp; j++)
330 for (
int d = 0; d < tar_ncomp; d++)
334 d*nsp*NE + i*nsp + j : i*nsp*
dim + d + j*
dim;
335 elem_dof_vals(j + d*nsp) = interp_vals(idx);
346 Vector elem_dof_vals(nsp*tar_ncomp);
348 for (
int i = 0; i < mesh_2.
GetNE(); i++)
352 for (
int j = 0; j < nsp; j++)
354 for (
int d = 0; d < tar_ncomp; d++)
358 d*nsp*NE + i*nsp + j : i*nsp*
dim + d + j*
dim;
359 elem_dof_vals(j*tar_ncomp+d) = interp_vals(idx);
378 cout <<
"Unable to connect to GLVis server at "
384 sout1 <<
"solution\n" << mesh_2 << func_target
385 <<
"window_title 'Target mesh and solution'"
386 <<
"window_geometry 600 0 600 600";
387 if (
dim == 2) { sout1 <<
"keys RmjAc"; }
388 if (
dim == 3) { sout1 <<
"keys mA\n"; }
394 ostringstream rho_name;
395 rho_name <<
"interpolated.gf";
396 ofstream rho_ofs(rho_name.str().c_str());
397 rho_ofs.precision(8);
398 func_target.
Save(rho_ofs);
405 if (func_source->
OwnFEC())
int Size() const
Return the logical size of the array.
Data type dense matrix using column-major storage.
void GetRow(int r, Vector &row) const
FindPointsGSLIB can robustly evaluate a GridFunction on an arbitrary collection of points....
void Setup(Mesh &m, const double bb_t=0.1, const double newt_tol=1.0e-12, const int npt_max=256)
virtual void Interpolate(const GridFunction &field_in, Vector &field_out)
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
virtual const char * Name() const
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
ElementTransformation * GetElementTransformation(int i) const
Returns ElementTransformation for the i-th element.
DofTransformation * GetElementVDofs(int i, Array< int > &vdofs) const
Returns indices of degrees of freedom for the i'th element. The returned indices are offsets into an ...
virtual const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th element in t...
Ordering::Type GetOrdering() const
Return the ordering method.
const FiniteElementCollection * FEColl() const
int GetVDim() const
Returns vector dimension.
Abstract class for all finite elements.
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...
const IntegrationRule & GetNodes() const
Get a const reference to the nodes of the element.
Geometry::Type GetGeomType() const
Returns the Geometry::Type of the reference element.
Class for grid function - Vector with associated FE space.
FiniteElementCollection * OwnFEC()
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
FiniteElementSpace * FESpace()
virtual void ProjectCoefficient(Coefficient &coeff)
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
Arbitrary order H1-conforming (continuous) finite elements.
Class for an integration rule - an Array of IntegrationPoint.
int GetNPoints() const
Returns the number of the points in the integration rule.
Arbitrary order "L2-conforming" discontinuous finite elements.
const FiniteElementSpace * GetNodalFESpace() const
int GetNE() const
Returns number of elements.
int Dimension() const
Dimension of the reference space used within the elements.
void GetNodes(Vector &node_coord) const
int GetNumGeometries(int dim) const
Return the number of geometries of the given dimension present in the mesh.
virtual void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1)
Set the curvature of the mesh nodes using the given polynomial degree.
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Arbitrary order H(curl)-conforming Nedelec finite elements.
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
void PrintUsage(std::ostream &out) const
Print the usage message.
void PrintOptions(std::ostream &out) const
Print the options.
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...
bool Good() const
Return true if the command line options were parsed successfully.
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
A general vector function coefficient.
void SetDataAndSize(real_t *d, int s)
Set the Vector data and size.
void SetSubVector(const Array< int > &dofs, const real_t value)
Set the entries listed in dofs to the given value.
int Size() const
Returns the size of the vector.
void SetSize(int s)
Resize the vector to size s.
real_t * GetData() const
Return a pointer to the beginning of the Vector data.
int open(const char hostname[], int port)
Open the socket stream on 'port' at 'hostname'.
double scalar_func(const Vector &x)
void vector_func(const Vector &p, Vector &F)
real_t p(const Vector &x, real_t t)