int main(int argc, char* argv[])
{
std::string func_name_x("x");
std::string func_name_y("y");
std::string func_name_z("z");
gsCmdLine cmd(
"Input is coordinate functions X(x,y,z), Y(x,y,z). Z(x,y,z), or"
"X(x,y), Y(x,y). Z(x,y), or X(x), Y(x). Z(x). This is controlled by"
"the dimension parameter d. The parameters take values in the interval [0,1]");
cmd.addInt("p", "degree", "this is the degree", p);
cmd.addInt("k", "knots", "This is the number of interior knots", k);
cmd.addInt("d", "dim", "this is the parametric dimension", d);
cmd.addInt("s", "samples", "this is samples for paraview", s);
cmd.addString("X", "f1", "The X-coordinate of the function", func_name_x);
cmd.addString("Y", "f2", "The Y-coordinate of the function", func_name_y);
cmd.addString("Z", "f3", "The Z-coordinate of the function", func_name_z);
try { cmd.getValues(argc,argv); } catch (int rv) { return rv; }
switch (d)
{
case 1:
break;
case 2:
break;
case 3:
break;
default:
{
gsWarn<<
"Dimension must be 1, 2 or 3.";
return 0;
}
};
gsInfo <<
"We are going to interpolate "<< func <<
"\n";
gsInfo <<
"Int. grid dim: "<< intGrid.dim() <<
"\n";
gsInfo <<
"Function values dim: "<< fValues.dim() <<
"\n";
gsInfo <<
"Result :"<< *interpolant <<
"\n";
mp.addPatch(
give(interpolant));
fd << mp;
std::string remark("Made by interpolation from function: ( ");
{
remark += func.expression(i);
remark += ( i==2 ? " )" : " , ");
}
fd.addComment( remark );
fd.save("interpolant_spline");
gsInfo <<
"Result saved as interpolant_spline.xml \n";
return 0;
}
A univariate B-spline basis.
Definition gsBSplineBasis.h:700
memory::unique_ptr< gsBasis > uPtr
Unique pointer for gsBasis.
Definition gsBasis.h:89
gsMatrix< T > anchors() const
Returns the anchor points that represent the members of the basis. There is exactly one anchor point ...
Definition gsBasis.h:437
virtual memory::unique_ptr< gsGeometry< T > > interpolateAtAnchors(gsMatrix< T > const &vals) const
Applies interpolation of values pts using the anchors as parameter points. May be reimplemented in de...
Definition gsBasis.hpp:267
Class for command-line argument parsing.
Definition gsCmdLine.h:57
This class represents an XML data tree which can be read from or written to a (file) stream.
Definition gsFileData.h:34
Class defining a multivariate (real or vector) function given by a string mathematical expression.
Definition gsFunctionExpr.h:52
memory::unique_ptr< gsGeometry > uPtr
Unique pointer for gsGeometry.
Definition gsGeometry.h:100
Class for representing a knot vector.
Definition gsKnotVector.h:80
A matrix with arbitrary coefficient type and fixed or dynamic size.
Definition gsMatrix.h:41
Container class for a set of geometry patches and their topology, that is, the interface connections ...
Definition gsMultiPatch.h:100
A tensor product B-spline basis.
Definition gsTensorBSplineBasis.h:37
Main header to be included by clients using the G+Smo library.
#define index_t
Definition gsConfig.h:32
#define gsWarn
Definition gsDebug.h:50
#define gsInfo
Definition gsDebug.h:43
unique_ptr< T > make_unique(T *x)
Definition gsMemory.h:198
The G+Smo namespace, containing all definitions for the library.
S give(S &x)
Definition gsMemory.h:266