12 #include "../mesh/mesh_headers.hpp"
60 MFEM_VERIFY(
T != NULL,
"invalid ElementTransformation");
68 double minDist = std::numeric_limits<double>::max();
72 for (
int i = 0; i < npts; ++i)
74 double dist = pt.
DistanceTo(physPts.GetColumn(i));
87 MFEM_VERIFY(
T != NULL,
"invalid ElementTransformation");
92 double minDist = std::numeric_limits<double>::max();
98 for (
int i = 0; i < npts; ++i)
105 double dist = dr.Norml2();
123 case 0: out <<
", ";
break;
124 case 1: out <<
"Newton: ";
break;
125 case 2: out <<
" ";
break;
130 case 0: out <<
"iter = " << std::setw(2) << int(val);
break;
131 case 1: out <<
"delta_ref = " << std::setw(11) << val;
break;
132 case 2: out <<
" err_phys = " << std::setw(11) << val;
break;
139 case 1: out <<
'\n';
break;
140 case 2: out <<
" (converged)\n";
break;
141 case 3: out <<
" (actual)\n";
break;
151 out << prefix <<
" = (";
152 for (
int j = 0; j < pt.
Size(); j++)
154 out << (j > 0 ?
", " :
"") << pt(j);
156 out <<
')' << suffix;
170 double xd[3], yd[3], dxd[3], dx_norm = -1.0, err_phys, real_dx_norm;
171 Vector x(xd, dim), y(yd, sdim), dx(dxd, dim);
172 bool hit_bdr =
false, prev_hit_bdr;
183 for (
int it = 0;
true; )
210 err_phys = y.Normlinf();
211 if (err_phys < phys_tol)
238 prev_xip.
Get(dxd, dim);
245 if (prev_hit_bdr && real_dx_norm <
ref_tol)
256 mfem::out <<
"Newton: *** stuck on boundary!\n";
273 prev_hit_bdr = hit_bdr;
285 default: MFEM_ABORT(
"invalid solver type");
317 mfem::out <<
"Newton: *** iteration did not converge!\n";
326 MFEM_VERIFY(
T != NULL,
"invalid ElementTransformation");
362 MFEM_ABORT(
"invalid initial guess type");
381 MFEM_ABORT(
"unknown Geometry::Type!");
384 int dof = FElem->
GetDof();
387 for (
int j = 0; j < dof; j++)
395 const DenseMatrix &IsoparametricTransformation::EvalJacobian()
398 "the IsoparametricTransformation has not been finalized;"
399 " call FinilizeTransformation() after setup");
404 if (dshape.
Width() > 0)
416 switch (FElem->
Space())
423 mfem_error(
"IsoparametricTransformation::OrderJ()");
430 switch (FElem->
Space())
437 mfem_error(
"IsoparametricTransformation::OrderW()");
452 return ((k-1)*(d-1)+(l-1));
454 return (k*(d-1)+(l-1));
457 mfem_error(
"IsoparametricTransformation::OrderGrad(...)");
467 FElem -> CalcShape(ip, shape);
468 PointMat.
Mult(shape, trans);
474 int dof, n,
dim, i, j, k;
483 for (j = 0; j < n; j++)
485 FElem -> CalcShape (ir.
IntPoint(j), shape);
486 for (i = 0; i <
dim; i++)
489 for (k = 0; k < dof; k++)
491 tr(i, j) += PointMat(i, k) * shape(k);
506 for (
int j = 0; j < matrix.
Width(); j++)
522 ip2.
Set(vec, v.Size());
531 for (i = 0; i < n; i++)
int GetNPoints() const
Returns the number of the points in the integration rule.
Abstract class for Finite Elements.
void Set(const double *p, const int dim)
void Get(double *p, const int dim) const
int GetDim() const
Returns the reference space dimension for the finite element.
Class for an integration rule - an Array of IntegrationPoint.
void SetSize(int s)
Resize the vector to size s.
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
RefinedGeometry * Refine(int Geom, int Times, int ETimes=1)
void CalcAdjugate(const DenseMatrix &a, DenseMatrix &adja)
int GetOrder() const
Returns the order of the finite element.
Data type dense matrix using column-major storage.
int Size() const
Returns the size of the vector.
int Space() const
Returns the type of space on each element.
double Normlinf() const
Returns the l_infinity norm of the vector.
const IntegrationPoint & GetCenter(int GeomType)
Return the center of the given Geometry::Type, GeomType.
PointFiniteElement PointFE
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
void SetType(const int t)
Set the Quadrature1D type of points to use for subdivision.
TriLinear3DFiniteElement HexahedronFE
const IntegrationRule & GetNodes() const
GeometryRefiner GlobGeometryRefiner
void GetColumn(int c, Vector &col) const
double DistanceTo(const double *p) const
Compute the Euclidean distance to another vector.
void CalcInverse(const DenseMatrix &a, DenseMatrix &inva)
static bool ProjectPoint(int GeomType, const IntegrationPoint &beg, IntegrationPoint &end)
Project a point end, onto the given Geometry::Type, GeomType.
int GetDof() const
Returns the number of degrees of freedom in the finite element.
void subtract(const Vector &x, const Vector &y, Vector &z)
void mfem_error(const char *msg)
void GetColumnReference(int c, Vector &col)
Linear3DFiniteElement TetrahedronFE
Linear2DFiniteElement TriangleFE
Class for integration point with weight.
int GetType() const
Get the Quadrature1D type of points used for subdivision.
void Mult(const double *x, double *y) const
Matrix vector multiplication.
static bool CheckPoint(int GeomType, const IntegrationPoint &ip)
Check if the given point is inside the given reference element.
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const =0
Evaluate the gradients of all shape functions of a scalar finite element in reference space at the gi...
BiLinear2DFiniteElement QuadrilateralFE
Linear1DFiniteElement SegmentFE
Tensor products of polynomials of order k.