|
IFEM
90A354
|
General utility classes and functions. More...
Classes | |
| class | Function |
| Base class for unary functions of arbitrary result and argument type. More... | |
| class | vector |
| A vector class with some added algebraic operations. More... | |
| class | matrixBase |
| Common base class for multi-dimensional (2D and 3D) matrices. More... | |
| class | matrix |
| Two-dimensional rectangular matrix with some algebraic operations. More... | |
| class | matrix3d |
| Three-dimensional rectangular matrix with some algebraic operations. More... | |
| class | matrix4d |
| Four-dimensional rectangular matrix with some algebraic operations. More... | |
| class | Function2 |
| Base class for binary function of arbitrary result and argument type. More... | |
| class | SpatialFunction |
| Base class for unary spatial function of arbitrary result type. More... | |
| class | LogStream |
| Logging stream class. More... | |
| class | Point |
| Class for representing points in 3D space. More... | |
| class | prof |
| Convenience class to profile the local scope. More... | |
| class | cmpInt |
| A helper class used to search for values in an integer map. More... | |
Functions | |
| template<class T > | |
| T | trunc (T v) |
| Truncate a value to zero when it is less than a given threshold. More... | |
| template<class T > | |
| std::istream & | operator>> (std::istream &s, vector< T > &X) |
| Read the vector X from the stream s. | |
| template<class T > | |
| std::ostream & | operator<< (std::ostream &s, const vector< T > &X) |
| Print the vector X to the stream s. | |
| template<class T > | |
| std::istream & | operator>> (std::istream &s, matrix< T > &A) |
| Read the matrix A from the stream s. | |
| template<class T > | |
| std::ostream & | operator<< (std::ostream &s, const matrix< T > &A) |
| Print the matrix A to the stream s. More... | |
| template<class T > | |
| void | writeMatlab (const char *label, const vector< T > &X, std::ostream &s=std::cout) |
| Print the vector X to the stream s in matlab format. | |
| template<class T > | |
| void | writeMatlab (const char *label, const matrix< T > &A, std::ostream &s=std::cout) |
| Print the matrix A to the stream s in matlab format. | |
| template<class T > | |
| std::ostream & | operator<< (std::ostream &s, const matrix3d< T > &A) |
| Print the 3D matrix A to the stream s. More... | |
| template<class T > | |
| std::ostream & | operator<< (std::ostream &s, const matrix4d< T > &A) |
| Print the 4D matrix A to the stream s. More... | |
| Vector | operator* (const Vector &X, Real c) |
| Multiplication of a vector and a scalar. More... | |
| Vector | operator* (Real c, const Vector &X) |
| Multiplication of a scalar and a vector. More... | |
| Vector | operator/ (const Vector &X, Real d) |
| Division of a vector by a scalar. More... | |
| Vector | operator+ (const Vector &X, const Vector &Y) |
| Addition of two vectors. More... | |
| Vector | operator- (const Vector &X, const Vector &Y) |
| Subtraction of two vectors. More... | |
| Matrix | operator* (const Matrix &A, Real c) |
| Multiplication of a matrix and a scalar. More... | |
| Matrix | operator* (Real c, const Matrix &A) |
| Multiplication of a scalar and a matrix. More... | |
| Real | operator* (const Vector &X, const Vector &Y) |
| Dot product of two vectors. More... | |
| RealArray | operator* (const Matrix &A, const Vector &X) |
| Multiplication of a matrix and a vector. More... | |
| RealArray | operator* (const Vector &X, const Matrix &A) |
| Multiplication of a vector and a matrix. More... | |
| Matrix | operator* (const Matrix &A, const Matrix &B) |
| Multiplication of two matrices. More... | |
| bool | transform (Matrix &A, const Matrix &Tn) |
| Congruence transformation of a symmetric matrix. More... | |
| bool | transform (Vector &V, const Matrix &Tn, bool transpose=false) |
| Congruence transformation of a vector. More... | |
| bool | invert (Matrix &A) |
| Inverts the square matrix A. | |
| bool | solve (Matrix &A, RealArray &b, std::vector< int > *iPivot=nullptr) |
| Solves the linear system of equations \( {\bf A x} = {\bf b} \). More... | |
| void | debugPrint (const char *label, const Vector &V, int level=2) |
| Debug print of given vector V, if SP_DEBUG ≥ level. | |
| Real | Jacobian (matrix< Real > &J, matrix< Real > &dNdX, const matrix< Real > &X, const matrix< Real > &dNdu, bool computeGradient=true) |
| Set up the Jacobian matrix of the coordinate mapping. More... | |
| Real | Jacobian (matrix< Real > &J, Vec3 &t, matrix< Real > &dNdX, const matrix< Real > &X, const matrix< Real > &dNdu, size_t tangent) |
| Set up the Jacobian matrix of the coordinate mapping along an edge. More... | |
| Real | Jacobian (matrix< Real > &J, Vec3 &n, matrix< Real > &dNdX, const matrix< Real > &X, const matrix< Real > &dNdu, size_t t1, size_t t2) |
| Set up the Jacobian matrix of the coordinate mapping on a boundary. More... | |
| bool | Hessian (matrix3d< Real > &H, matrix3d< Real > &d2NdX2, const matrix< Real > &Ji, const matrix< Real > &X, const matrix3d< Real > &d2Ndu2, const matrix< Real > &dNdX, bool geoMapping=true) |
| Set up the Hessian matrix of the coordinate mapping. More... | |
| void | Hessian (const matrix3d< Real > &Hess, matrix< Real > &H) |
| Convert a Hessian from a matrix3d to a matrix assuming symmetry. | |
| void | getGmat (const matrix< Real > &Ji, const Real *du, matrix< Real > &G) |
| Compute the stabilization matrix G from the Jacobian inverse. More... | |
| bool | Hessian2 (matrix4d< Real > &d3NdX3, const matrix< Real > &Ji, const matrix4d< Real > &d3Ndu3) |
| Set up third-order derivatives of the coordinate mapping. More... | |
| void | JacobianGradient (const matrix< Real > &dudX, const matrix3d< Real > &d2Xdu2, std::vector< matrix< Real >> &dJdX) |
| Calculates the derivatives of the Jacobian of the coordinate mapping. More... | |
| void | detJacGradient (const matrix< Real > &J, const matrix< Real > &Ji, const matrix3d< Real > &H, std::vector< Real > &ddet) |
| Calculates the derivatives of determinant of the Jacobian. More... | |
| const RealFunc * | parseRealFunc (char *cline, Real A=Real(1), bool print=true) |
| Creates a scalar-valued function by parsing a character string. More... | |
| ScalarFunc * | parseTimeFunc (const char *func, const std::string &type="expression", Real eps=Real(1.0e-8)) |
| Creates a time function by parsing a character string. More... | |
| VecTimeFunc * | parseVecTimeFunc (const char *func, const std::string &type) |
| Creates a vector-valued time function by parsing a char string. More... | |
| RealFunc * | parseRealFunc (const std::string &func, const std::string &type="expression", bool print=true) |
| Creates a scalar-valued function by parsing a character string. More... | |
| RealFunc * | parseExprRealFunc (const std::string &function, bool autodiff) |
| Creates a scalar-valued function by parsing a character string. More... | |
| VecFunc * | parseExprVecFunc (const std::string &function, bool autodiff) |
| Creates a Vector-valued function by parsing a character string. More... | |
| IntFunc * | parseIntFunc (const std::string &func, const std::string &type="expression") |
| Creates a scalar-valued int function by parsing a character string. More... | |
| VecFunc * | parseVecFunc (const std::string &func, const std::string &type="expression", const std::string &variables="") |
| Creates a vector-valued function by parsing a character string. More... | |
| TensorFunc * | parseTensorFunc (const std::string &func, const std::string &type) |
| Creates a tensor-valued function by parsing a character string. More... | |
| TractionFunc * | parseTracFunc (const std::string &func, const std::string &type="expression", int dir=0) |
| Creates a vector-valued function defining a surface traction. More... | |
| TractionFunc * | parseTracFunc (const tinyxml2::XMLElement *elem) |
| Creates a vector-valued function defining a surface traction. More... | |
| size_t | Pascal (int p, unsigned short int nsd) |
| Returns the number of monomials in Pascal's triangle. More... | |
| void | Pascal (int p, Real x, Real y, std::vector< Real > &phi) |
| Evaluates the monomials of Pascal's triangle in 2D for order p. | |
| void | Pascal (int p, Real x, Real y, Real z, std::vector< Real > &phi) |
| Evaluates the monomials of Pascal's triangle in 3D for order p. | |
| Real | Pos (Real x) |
| Evaluates the positive ramp function \((x+|x|)/2\). | |
| Real | Neg (Real x) |
| Evaluates the negative ramp function \((x-|x|)/2\). | |
| Real | Rad (Real x) |
| Converts from degrees to radians. | |
| bool | parseIntegers (std::vector< int > &values, const char *argv) |
| Parses a character string into an integer or an integer range. More... | |
| bool | parseVec (Vec3 &value, const char *argv) |
| Parses a character string into a point vector. More... | |
| bool | parseKnots (std::vector< Real > &xi) |
| Parses a (possibly graded) sequence of knot values. More... | |
| char * | readLine (std::istream &is) |
| Reads one line, ignoring comment lines and leading blanks. More... | |
| bool | ignoreComments (std::istream &is) |
| Ignores comment lines and blank lines from an input stream. More... | |
| int | getAttribute (const tinyxml2::XMLElement *xml, const char *att, bool &val) |
| Extracts a boolean attribute value from the specified XML-element. More... | |
| int | getAttribute (const tinyxml2::XMLElement *xml, const char *att, int &val) |
| Extracts an integer attribute value from the specified XML-element. More... | |
| int | getAttribute (const tinyxml2::XMLElement *xml, const char *att, char &val, bool useIntValue=true) |
| Extracts a char attribute value from the specified XML-element. More... | |
| int | getAttribute (const tinyxml2::XMLElement *xml, const char *att, size_t &val) |
| Extracts a size_t attribute value from the specified XML-element. More... | |
| int | getAttribute (const tinyxml2::XMLElement *xml, const char *att, Real &val) |
| Extracts a real attribute value from the specified XML-element. More... | |
| int | getAttribute (const tinyxml2::XMLElement *xml, const char *att, Vec3 &val, int ncomp=0) |
| Extracts a vector attribute value from the specified XML-element. More... | |
| int | getAttribute (const tinyxml2::XMLElement *xml, const char *att, std::vector< int > &val) |
| Extracts an integer vector attribute from specified XML-element. More... | |
| int | getAttribute (const tinyxml2::XMLElement *xml, const char *att, std::string &val, bool toLower=false) |
| Extracts a string attribute value from the specified XML-element. More... | |
| const char * | getValue (const tinyxml2::XMLNode *xml, const char *tag) |
| Returns the value (if any) of the specified XML-node. More... | |
| bool | parseKnots (const tinyxml2::XMLNode *xml, std::vector< Real > &xi) |
| Parses a sequence of knot values from the specified XML-node. More... | |
| bool | renumber (int &num, int &runner, IntMap &old2new) |
| Transforms the integer value num into a unique range. More... | |
| bool | renumber (int &num, const IntMap &old2new, bool msg=false) |
| Transforms the integer value num according to a given mapping. More... | |
| int | gather (const std::vector< int > &index, size_t nr, const std::vector< Real > &in, std::vector< Real > &out, size_t offset_in=0, size_t nskip=0) |
| Compresses the rows of a 2D array based on given scatter indices. More... | |
| int | gather (const std::vector< int > &index, size_t nr, const std::vector< Real > &in, utl::matrix< Real > &out, size_t offset_in=0, size_t nskip=0) |
| Compresses the rows of a 2D array based on given scatter indices. More... | |
| int | gather (const std::vector< int > &index, size_t ir, size_t nr, const std::vector< Real > &in, std::vector< Real > &out, size_t offset_in=0, int shift_idx=0, size_t nskip=0) |
| Compresses a row of a 2D array based on given scatter indices. More... | |
| void | interleave (const std::vector< Real > &v1, const std::vector< Real > &v2, std::vector< Real > &out, size_t n1=1, size_t n2=1) |
| Merges two arrays of nodal values into one array by interleaving. More... | |
| size_t | find_closest (const std::vector< Real > &a, Real v) |
| Searches for a real value in an ordered array of reals. More... | |
| std::string | adjustRight (size_t width, const std::string &s, const std::string &suffix=" : ") |
| Right-justifies the input string to the given total width. | |
| std::set< int > | getDigits (int value) |
| Splits an integer into its (unique) digits in ascending order. | |
| int | getDirs (int ncmp) |
| Makes a direction integer for given number of components. | |
| IntMap::const_iterator | findValue (const IntMap &iMap, int iVal) |
| Returns a const iterator to the entry with value iVal. | |
| int | findKey (const IntMap &iMap, int iVal) |
| Returns the key corresponding to the value iVal. More... | |
| int | findIndex (const std::vector< int > &iVec, int iVal) |
| Returns the vector index corresponding to the value iVal. More... | |
| void | merge (std::vector< int > &a1, const std::vector< int > &a2) |
| Merges integer array a2 into array a1. More... | |
| void | merge (std::vector< Real > &a1, const std::vector< Real > &a2, const std::vector< int > &k1, const std::vector< int > &k2) |
| Merges real array a2 into array a1 based on array indices. More... | |
| Real | round (Real value, int pr) |
| Rounds off value down to pr significant digits. | |
Variables | |
| const char | RETAIN = 2 |
| Flag for vector::resize() method telling it to retain its content. | |
| double | zero_print_tol = 1.0e-6 |
| Zero tolerance for printing numbers. | |
| int | nval_per_line = 6 |
| Number of values to print per line. | |
| Profiler * | profiler = nullptr |
| Pointer to the one and only profiler object. | |
General utility classes and functions.
| void utl::detJacGradient | ( | const matrix< Real > & | J, |
| const matrix< Real > & | Ji, | ||
| const matrix3d< Real > & | H, | ||
| std::vector< Real > & | ddet | ||
| ) |
Calculates the derivatives of determinant of the Jacobian.
| [in] | J | Jacobian of the geometry mapping |
| [in] | Ji | Inverse jacobian of the geometry mapping |
| [in] | H | Hessian of the geometry mapping wrt parameters |
| [out] | ddet | Derivatives of the determinant |
References utl::matrix< T >::multiply(), Real, and utl::matrix< T >::rows().
Referenced by MxFiniteElement::piolaGradient().
Searches for a real value in an ordered array of reals.
| [in] | a | The array of real values |
| [in] | v | The value to search for |
Referenced by ASMs1DLag::evalPoint(), ASMs2DLag::evalPoint(), and ASMs3DLag::evalPoint().
| int utl::findIndex | ( | const std::vector< int > & | iVec, |
| int | iVal | ||
| ) |
Returns the vector index corresponding to the value iVal.
If not in the vector, -1 is returned.
Referenced by ASMbase::addLagrangeMultipliers(), ASMu2D::constrainEdge(), ASMbase::getAge(), ASMbase::getElmIndex(), ASMLRSpline::getFunctionsForElements(), ASMs2D::getNodeIndex(), ASMs3D::getNodeIndex(), ASMbase::getNodeIndex(), ASMu2DIB::integrate(), ASMbase::isElementInPartition(), ASMu1DLag::isInElementSet(), ASMu2DLag::isInElementSet(), ASMs1D::isInNodeSet(), ASMu2DLag::isInNodeSet(), ASMu2DLag::parseElemBox(), and ASMu2DLag::parseNodeBox().
| int utl::findKey | ( | const IntMap & | iMap, |
| int | iVal | ||
| ) |
Returns the key corresponding to the value iVal.
If not in the map, the value iVal is returned.
References findValue().
Referenced by ASMbase::updateDirichlet().
| int utl::gather | ( | const std::vector< int > & | index, |
| size_t | ir, | ||
| size_t | nr, | ||
| const std::vector< Real > & | in, | ||
| std::vector< Real > & | out, | ||
| size_t | offset_in = 0, |
||
| int | shift_idx = 0, |
||
| size_t | nskip = 0 |
||
| ) |
Compresses a row of a 2D array based on given scatter indices.
| [in] | index | Scatter indices of the columns that should be retained |
| [in] | ir | Index of the row to compress |
| [in] | nr | Number of rows in the 2D array |
| [in] | in | The input array stored column-wise in a 1D array |
| [out] | out | The output array stored column-wise in a 1D array |
| [in] | offset_in | Optional start offset for the in vector |
| [in] | shift_idx | Optional constant shift in the scatter indices |
| [in] | nskip | If nonzero, skip the last nskip entries in index |
| int utl::gather | ( | const std::vector< int > & | index, |
| size_t | nr, | ||
| const std::vector< Real > & | in, | ||
| std::vector< Real > & | out, | ||
| size_t | offset_in = 0, |
||
| size_t | nskip = 0 |
||
| ) |
Compresses the rows of a 2D array based on given scatter indices.
| [in] | index | Scatter indices of the columns that should be retained |
| [in] | nr | Number of rows in the 2D array |
| [in] | in | The input array stored column-wise in a 1D array |
| [out] | out | The output array stored column-wise in a 1D array |
| [in] | offset_in | Optional start offset for the in vector |
| [in] | nskip | If nonzero, skip the last nskip entries in index |
References Real.
Referenced by IntegrandBase::evalSol1(), ASMs2Dmx::evalSolution(), ASMs3Dmx::evalSolution(), ASMs2D::evalSolution(), ASMs3D::evalSolution(), ASMu2D::evalSolution(), ASMu3D::evalSolution(), ASMs1D::evalSolution(), ASMs2Dmx::evalSolutionPiola(), SplineFields1D::gradFE(), SplineFields2D::gradFE(), SplineFields2Dmx::gradFE(), SplineFields3D::gradFE(), SplineFields3Dmx::gradFE(), SplineField2D::gradFE(), SplineField3D::gradFE(), SplineField2D::hessianFE(), SplineField3D::hessianFE(), SplineFields1D::hessianFE(), SplineFields2D::hessianFE(), SplineFields2Dmx::hessianFE(), SplineFields3D::hessianFE(), SplineFields3Dmx::hessianFE(), IntegrandBase::initElement1(), IntegrandBase::initElement2(), NormBase::initProjection(), SplineField2D::valueFE(), SplineField3D::valueFE(), SplineFields1D::valueFE(), SplineFields2D::valueFE(), SplineFields2Dmx::valueFE(), SplineFields3D::valueFE(), SplineFields3Dmx::valueFE(), SplineField2D::valueGrid(), and SplineField3D::valueGrid().
| int utl::gather | ( | const std::vector< int > & | index, |
| size_t | nr, | ||
| const std::vector< Real > & | in, | ||
| utl::matrix< Real > & | out, | ||
| size_t | offset_in = 0, |
||
| size_t | nskip = 0 |
||
| ) |
Compresses the rows of a 2D array based on given scatter indices.
| [in] | index | Scatter indices of the columns that should be retained |
| [in] | nr | Number of rows in the 2D array |
| [in] | in | The input array stored column-wise in a 1D array |
| [out] | out | The output array stored as a 2D matrix |
| [in] | offset_in | Optional start offset for the in vector |
| [in] | nskip | If nonzero, skip the last nskip entries in index |
References utl::matrix< T >::fillColumn(), Real, and utl::matrix< T >::resize().
| int utl::getAttribute | ( | const tinyxml2::XMLElement * | xml, |
| const char * | att, | ||
| bool & | val | ||
| ) |
Extracts a boolean attribute value from the specified XML-element.
| [in] | xml | Pointer to XML-element to extract from |
| [in] | att | The attribute tag |
| [out] | val | The attribute value |
Referenced by AnaSol::AnaSol(), MultiPatchModelGenerator1D::createG2(), MultiPatchModelGenerator2D::createG2(), DefaultGeometry1D::createG2(), DefaultGeometry2D::createG2(), MultiPatchModelGenerator3D::createG2(), DefaultGeometry3D::createG2(), MultiPatchModelGenerator1D::createGeometry(), MultiPatchModelGenerator2D::createGeometry(), MultiPatchModelGenerator3D::createGeometry(), ModelGenerator::createGeometry(), MultiPatchModelGenerator1D::MultiPatchModelGenerator1D(), MultiPatchModelGenerator2D::MultiPatchModelGenerator2D(), MultiPatchModelGenerator3D::MultiPatchModelGenerator3D(), DataExporter::OnControl(), SIMNodalConstraint< Dim >::parse(), SIMSolverTS< T1 >::parse(), HasGravityBase::parse(), AdaptiveSetup::parse(), EigenModeSIM::parse(), NewmarkSIM::parse(), NonLinSIM::parse(), SIM1D::parse(), SIM2D::parse(), SIMinput::parse(), SIMsupel::parse(), TimeStep::parse(), TextureProperties::parse(), SIMSemi3D< PlaneSolver >::parse(), SIMargsBase::parse(), SIM1D::parseBCTag(), SIM2D::parseBCTag(), SIM3D::parseBCTag(), SIMinput::parseBCTag(), SIMmultiCpl::parseConnection(), SIMoptions::parseConsoleTag(), SIMoptions::parseDiscretizationTag(), SIMinput::parseDualTag(), AnaSol::parseExpressionFunctions(), AnaSol::parseFieldFunctions(), SIM1D::parseGeometryTag(), SIM2D::parseGeometryTag(), SIM3D::parseGeometryTag(), SIMinput::parseICTag(), SIMCoupledSI< T1, T2 >::parseIterations(), parseKnots(), SIMinput::parseLinSolTag(), SIMinput::parseMaterialSet(), SIMinput::parseOutputTag(), SIMoptions::parseOutputTag(), SIMoutput::parseOutputTag(), SIMmodal::parseParams(), SIMinput::parsePatchList(), SIMinput::parsePeriodic(), SIMoptions::parseRestartTag(), SIMinput::parseTopologySet(), parseTracFunc(), LinSolParams::BlockParams::read(), ASM::readXML(), and ModelGenerator::topologySets().
| int utl::getAttribute | ( | const tinyxml2::XMLElement * | xml, |
| const char * | att, | ||
| char & | val, | ||
| bool | useIntValue = true |
||
| ) |
Extracts a char attribute value from the specified XML-element.
| [in] | xml | Pointer to XML-element to extract from |
| [in] | att | The attribute tag |
| [out] | val | The attribute value |
| [in] | useIntValue | If true, convert the value to an integer |
| int utl::getAttribute | ( | const tinyxml2::XMLElement * | xml, |
| const char * | att, | ||
| int & | val | ||
| ) |
Extracts an integer attribute value from the specified XML-element.
| [in] | xml | Pointer to XML-element to extract from |
| [in] | att | The attribute tag |
| [out] | val | The attribute value |
| int utl::getAttribute | ( | const tinyxml2::XMLElement * | xml, |
| const char * | att, | ||
| Real & | val | ||
| ) |
Extracts a real attribute value from the specified XML-element.
| [in] | xml | Pointer to XML-element to extract from |
| [in] | att | The attribute tag |
| [out] | val | The attribute value |
| int utl::getAttribute | ( | const tinyxml2::XMLElement * | xml, |
| const char * | att, | ||
| size_t & | val | ||
| ) |
Extracts a size_t attribute value from the specified XML-element.
| [in] | xml | Pointer to XML-element to extract from |
| [in] | att | The attribute tag |
| [out] | val | The attribute value |
| int utl::getAttribute | ( | const tinyxml2::XMLElement * | xml, |
| const char * | att, | ||
| std::string & | val, | ||
| bool | toLower = false |
||
| ) |
Extracts a string attribute value from the specified XML-element.
| [in] | xml | Pointer to XML-element to extract from |
| [in] | att | The attribute tag |
| [out] | val | The attribute value |
| [in] | toLower | If true, convert return string to lower case |
| int utl::getAttribute | ( | const tinyxml2::XMLElement * | xml, |
| const char * | att, | ||
| std::vector< int > & | val | ||
| ) |
Extracts an integer vector attribute from specified XML-element.
| [in] | xml | Pointer to XML-element to extract from |
| [in] | att | The attribute tag |
| [out] | val | The attribute value |
References parseIntegers().
| int utl::getAttribute | ( | const tinyxml2::XMLElement * | xml, |
| const char * | att, | ||
| Vec3 & | val, | ||
| int | ncomp = 0 |
||
| ) |
Extracts a vector attribute value from the specified XML-element.
| [in] | xml | Pointer to XML-element to extract from |
| [in] | att | The attribute tag |
| [out] | val | The attribute value |
| [in] | ncomp | Maximum number of components to read |
If fewer than ncomp components specified, use the last value for the rest. If ncomp is zero, use the value zero for the missing components, if any.
Compute the stabilization matrix G from the Jacobian inverse.
| [in] | Ji | The inverse of the Jacobian matrix |
| [in] | du | Element lengths in each parametric direction |
| [out] | G | The stabilization matrix (used in CFD simulators) |
References utl::matrix< T >::cols(), Real, and utl::matrix< T >::resize().
Referenced by ASMs2D::integrate(), ASMs2Dmx::integrate(), ASMs3D::integrate(), ASMs3Dmx::integrate(), ASMu2D::integrate(), ASMu2Dmx::integrate(), ASMu3D::integrate(), and ASMu3Dmx::integrate().
| const char * utl::getValue | ( | const tinyxml2::XMLNode * | xml, |
| const char * | tag | ||
| ) |
Returns the value (if any) of the specified XML-node.
| [in] | xml | Pointer to XML-node to extract the value from |
| [in] | tag | The name of the XML-element to extract the value from |
This method accepts two alternative ways of specifying the value myValue :
<name>myValue</name>
and
<name value="myValue"/>
Referenced by SIMSolverTS< T1 >::parse(), AdaptiveSetup::parse(), EigenModeSIM::parse(), NewmarkSIM::parse(), NonLinSIM::parse(), SIMadmin::parse(), SIMoptions::parseEigSolTag(), SIMoptions::parseOutputTag(), SIMoutput::parseOutputTag(), LinSolParams::BlockParams::read(), and LinSolParams::read().
| bool utl::Hessian | ( | matrix3d< Real > & | H, |
| matrix3d< Real > & | d2NdX2, | ||
| const matrix< Real > & | Ji, | ||
| const matrix< Real > & | X, | ||
| const matrix3d< Real > & | d2Ndu2, | ||
| const matrix< Real > & | dNdX, | ||
| bool | geoMapping = true |
||
| ) |
Set up the Hessian matrix of the coordinate mapping.
| [out] | H | The Hessian matrix |
| [out] | d2NdX2 | Second order derivatives of basis functions, w.r.t. X |
| [in] | Ji | The inverse of the Jacobian matrix |
| [in] | X | Matrix of element nodal coordinates |
| [in] | d2Ndu2 | Second order derivatives of basis functions |
| [in] | dNdX | First order derivatives of basis functions |
| [in] | geoMapping | If true, calculate geometry mapping |
If geoMapping is true, H is output, else H is input and assumed to be already calculated in a previous call.
References utl::matrixBase< T >::clear(), utl::matrix< T >::cols(), utl::matrixBase< T >::empty(), utl::matrix3d< T >::multiply(), PROFILE4, Real, utl::matrix3d< T >::resize(), and utl::matrix< T >::rows().
Referenced by SplineField::evalBasis(), LRSplineField::evalBasis(), SplineField::evalMapping(), LRSplineField::evalMapping(), ASMs2D::evalSolution(), ASMs3D::evalSolution(), ASMu2D::evalSolution(), ASMu3D::evalSolution(), ASMs1D::evalSolution(), MxFiniteElement::Hessian(), SplineFields1D::hessianFE(), ASMs1D::integrate(), ASMs2D::integrate(), ASMs3D::integrate(), ASMu2D::integrate(), and ASMu3D::integrate().
| bool utl::Hessian2 | ( | matrix4d< Real > & | d3NdX3, |
| const matrix< Real > & | Ji, | ||
| const matrix4d< Real > & | d3Ndu3 | ||
| ) |
Set up third-order derivatives of the coordinate mapping.
| [out] | d3NdX3 | Third order derivatives of basis functions, w.r.t. X |
| [in] | Ji | The inverse of the Jacobian matrix |
| [in] | d3Ndu3 | Third order derivatives of basis functions |
References utl::matrix< T >::cols(), utl::matrixBase< T >::dim(), PROFILE4, utl::matrix4d< T >::resize(), and utl::matrix< T >::rows().
Referenced by ASMs2D::evalSolution(), ASMu2D::evalSolution(), ASMs1D::evalSolution(), ASMs1D::integrate(), ASMs2D::integrate(), and ASMu2D::integrate().
| bool utl::ignoreComments | ( | std::istream & | is | ) |
Ignores comment lines and blank lines from an input stream.
| is | File stream to read from |
| void utl::interleave | ( | const std::vector< Real > & | v1, |
| const std::vector< Real > & | v2, | ||
| std::vector< Real > & | out, | ||
| size_t | n1 = 1, |
||
| size_t | n2 = 1 |
||
| ) |
Merges two arrays of nodal values into one array by interleaving.
| [in] | v1 | The first array |
| [in] | v2 | The second array |
| [out] | out | The output array |
| [in] | n1 | Number of entries per node in first array |
| [in] | n2 | Number of entries per node in second array |
Referenced by ASMu2D::regularInterpolation(), and ASMu3D::regularInterpolation().
| Real utl::Jacobian | ( | matrix< Real > & | J, |
| matrix< Real > & | dNdX, | ||
| const matrix< Real > & | X, | ||
| const matrix< Real > & | dNdu, | ||
| bool | computeGradient = true |
||
| ) |
Set up the Jacobian matrix of the coordinate mapping.
| [out] | J | The inverse of the Jacobian matrix |
| [out] | dNdX | First order derivatives of basis functions, w.r.t. X |
| [in] | X | Matrix of element nodal coordinates |
| [in] | dNdu | First order derivatives of basis functions |
| [in] | computeGradient | If false, skip calculation of dNdX |
References utl::matrixBase< T >::clear(), utl::matrix< T >::cols(), epsZ, utl::matrix< T >::getColumn(), utl::matrix< T >::inverse(), Vec3::length(), utl::matrix< T >::multiply(), Real, and utl::matrix< T >::rows().
Referenced by ASMs1D::assembleL2matrices(), ASMs2DLag::assembleL2matrices(), ASMs3DLag::assembleL2matrices(), ASMu2D::assembleL2matrices(), ASMu3D::assembleL2matrices(), SplineField::evalMapping(), LRSplineField::evalMapping(), ASMs2DmxLag::evalSolution(), ASMs3DmxLag::evalSolution(), ASMs2DLag::evalSolution(), ASMs3DLag::evalSolution(), ASMs2DTri::evalSolution(), ASMs1DLag::evalSolution(), ASMs1DSpec::evalSolution(), ASMs2DSpec::evalSolution(), ASMs3DSpec::evalSolution(), ASMs2D::evalSolution(), ASMs3D::evalSolution(), ASMu2D::evalSolution(), ASMu3D::evalSolution(), ASMs1D::evalSolution(), ASMu2D::evaluateBasisNurbs(), LagrangeFields2D::gradFE(), LagrangeFields3D::gradFE(), SplineFields1D::gradFE(), LagrangeField3D::gradFE(), LagrangeField2D::gradFE(), SplineFields1D::hessianFE(), ASMs1D::integrate(), ASMs1DLag::integrate(), ASMs1DSpec::integrate(), ASMs2D::integrate(), ASMs2DSpec::integrate(), ASMs2DTri::integrate(), ASMs3D::integrate(), ASMs3DLag::integrate(), ASMs3DSpec::integrate(), ASMu2D::integrate(), ASMu3D::integrate(), ASMu3Dmx::integrate(), ASMs2DLag::integrate(), ASMs3D::integrateEdge(), ASMs3DLag::integrateEdge(), ASMs3DSpec::integrateEdge(), ASMs2DLag::integrateElm(), and MxFiniteElement::Jacobian().
| Real utl::Jacobian | ( | matrix< Real > & | J, |
| Vec3 & | n, | ||
| matrix< Real > & | dNdX, | ||
| const matrix< Real > & | X, | ||
| const matrix< Real > & | dNdu, | ||
| size_t | t1, | ||
| size_t | t2 | ||
| ) |
Set up the Jacobian matrix of the coordinate mapping on a boundary.
| [out] | J | The inverse of the Jacobian matrix |
| [out] | n | Outward-directed unit normal vector on the boundary |
| [out] | dNdX | 1st order derivatives of basis functions, w.r.t. X |
| [in] | X | Matrix of element nodal coordinates |
| [in] | dNdu | First order derivatives of basis functions |
| [in] | t1 | First parametric tangent direction of the boundary |
| [in] | t2 | Second parametric tangent direction of the boundary |
References utl::matrixBase< T >::clear(), utl::matrix< T >::cols(), Vec3::cross(), epsZ, utl::matrix< T >::getColumn(), utl::matrix< T >::inverse(), utl::matrix< T >::multiply(), Vec3::normalize(), Real, and utl::matrix< T >::rows().
| Real utl::Jacobian | ( | matrix< Real > & | J, |
| Vec3 & | t, | ||
| matrix< Real > & | dNdX, | ||
| const matrix< Real > & | X, | ||
| const matrix< Real > & | dNdu, | ||
| size_t | tangent | ||
| ) |
Set up the Jacobian matrix of the coordinate mapping along an edge.
| [out] | J | The inverse of the Jacobian matrix |
| [out] | t | Unit tangent vector along the edge |
| [out] | dNdX | 1st order derivatives of basis functions, w.r.t. X |
| [in] | X | Matrix of element nodal coordinates |
| [in] | dNdu | First order derivatives of basis functions |
| [in] | tangent | Parametric tangent direction along the edge |
References utl::matrixBase< T >::clear(), epsZ, utl::matrix< T >::getColumn(), utl::matrix< T >::inverse(), utl::matrix< T >::multiply(), Vec3::normalize(), and Real.
| void utl::JacobianGradient | ( | const matrix< Real > & | dudX, |
| const matrix3d< Real > & | d2Xdu2, | ||
| std::vector< matrix< Real >> & | dJdX | ||
| ) |
Calculates the derivatives of the Jacobian of the coordinate mapping.
| [in] | dudX | Derivatives of the geometry basis |
| [in] | d2Xdu2 | Second order derivatives of the geometry basis |
| [out] | dJdX | Derivatives of the Jacobian wrt physical coordinates |
References utl::matrix< T >::add(), utl::matrix< T >::resize(), and utl::matrix< T >::rows().
Referenced by MxFiniteElement::piolaGradient().
| void utl::merge | ( | std::vector< int > & | a1, |
| const std::vector< int > & | a2 | ||
| ) |
Merges integer array a2 into array a1.
Does not require the arrays to be sorted. The values of a2 not already in a1 are appended to a1.
Referenced by ASMs2D::addInterfaceElms(), and ASMs2D::integrate().
| void utl::merge | ( | std::vector< Real > & | a1, |
| const std::vector< Real > & | a2, | ||
| const std::vector< int > & | k1, | ||
| const std::vector< int > & | k2 | ||
| ) |
Merges real array a2 into array a1 based on array indices.
Does not require the arrays to be sorted. The values of a2 not already in a1 are appended to a1.
Multiplication of two matrices.
References utl::matrix< T >::multiply().
Multiplication of a matrix and a vector.
References utl::matrix< T >::multiply().
Multiplication of a matrix and a scalar.
References utl::matrix< T >::multiply().
Multiplication of a vector and a matrix.
References utl::matrix< T >::multiply().
Multiplication of a vector and a scalar.
Multiplication of a scalar and a matrix.
Multiplication of a scalar and a vector.
Addition of two vectors.
References utl::vector< T >::add().
Subtraction of two vectors.
References utl::vector< T >::add(), and Real.
| std::ostream& utl::operator<< | ( | std::ostream & | s, |
| const matrix3d< T > & | A | ||
| ) |
Print the 3D matrix A to the stream s.
The matrix is printed as a set of 2D sub-matrices based on the first two indices.
References utl::matrixBase< T >::dim(), utl::matrixBase< T >::empty(), utl::matrix< T >::fill(), utl::matrixBase< T >::ptr(), and utl::matrixBase< T >::size().
| std::ostream& utl::operator<< | ( | std::ostream & | s, |
| const matrix4d< T > & | A | ||
| ) |
Print the 4D matrix A to the stream s.
The matrix is printed as a set of 2D sub-matrices based on the first two indices.
References utl::matrixBase< T >::dim(), utl::matrixBase< T >::empty(), utl::matrix< T >::fill(), utl::matrixBase< T >::ptr(), and utl::matrixBase< T >::size().
| std::ostream& utl::operator<< | ( | std::ostream & | s, |
| const matrix< T > & | A | ||
| ) |
Print the matrix A to the stream s.
If the matrix is symmetric, only the upper triangular part of the matrix is printed, with the diagonal elements in the first column. The global variable zero_print_tol is used as a tolerance when checking whether the matrix is symmetric or not.
References utl::matrix< T >::cols(), utl::matrix< T >::isSymmetric(), utl::matrix< T >::rows(), trunc(), and zero_print_tol.
| RealFunc * utl::parseExprRealFunc | ( | const std::string & | function, |
| bool | autodiff | ||
| ) |
Creates a scalar-valued function by parsing a character string.
| [in] | function | Function expression |
| [in] | autodiff | If true, auto-differentiation is enabled |
The implementation of this method is in the file ExprFunctions.C for encapsulation of the autodiff package.
| VecFunc * utl::parseExprVecFunc | ( | const std::string & | function, |
| bool | autodiff | ||
| ) |
Creates a Vector-valued function by parsing a character string.
| [in] | function | Function expression |
| [in] | autodiff | If true, auto-differentiation is enabled |
The implementation of this method is in the file ExprFunctions.C for encapsulation of the autodiff package.
| bool utl::parseIntegers | ( | std::vector< int > & | values, |
| const char * | argv | ||
| ) |
Parses a character string into an integer or an integer range.
| values | The integer value(s) is/are appended to this vector | |
| [in] | argv | Character string with integer data |
An integer range is recognised through the syntax i:j.
Referenced by getAttribute(), EigenModeSIM::parse(), ASMu1DLag::parseElemSet(), ASMu2DLag::parseElemSet(), ASMs1D::parseNodeSet(), ASMsupel::parseNodeSet(), and ASMu2DLag::parseNodeSet().
| IntFunc * utl::parseIntFunc | ( | const std::string & | func, |
| const std::string & | type = "expression" |
||
| ) |
Creates a scalar-valued int function by parsing a character string.
| [in] | func | Character string to parse function definition from |
| [in] | type | Function definition type flag |
References IFEM::cout, utl::Function< Arg, Result >::evaluate(), and Real.
| bool utl::parseKnots | ( | const tinyxml2::XMLNode * | xml, |
| std::vector< Real > & | xi | ||
| ) |
Parses a sequence of knot values from the specified XML-node.
| [in] | xml | Pointer to XML-node to extract from |
| xi | The knot value(s) is/are appended to this vector |
References parseKnots().
| bool utl::parseKnots | ( | std::vector< Real > & | xi | ) |
Parses a (possibly graded) sequence of knot values.
| xi | The knot value(s) is/are appended to this vector |
The method uses the strtok function to parse a sequence of space-separated numbers and assumes that (at least one) call to strtok with the first argument different from nullptr has been made before invoking this method.
If the first character of the first token is a digit, it is assumed that the knot values are listed explicitly by the remaining tokens. Otherwise, various mesh grading schemes are assumed with the following syntax:
In each case the starting and ending knot values, xi1 and xi2, are optional. If not specified, 0.0 and 1.0 is assumed. The 4th scheme was proposed by Tymofiy Gerasimov 14.12.2016, and is described in the file mesh-grading-function.pdf located in the doc folder.
Referenced by SIM1D::parse(), SIM2D::parse(), SIM3D::parse(), SIM1D::parseGeometryTag(), SIM2D::parseGeometryTag(), SIM3D::parseGeometryTag(), and parseKnots().
Creates a scalar-valued function by parsing a character string.
The functions are assumed on the general form
\[ f({\bf X},t) = A * g({\bf X}) * h(t) \]
The character string cline is assumed to contain first the definition of the spatial function g( X ) and then the time function h(t). Either of the two components may be omitted, for creating a space-function constant in time, or a time function constant in space.
References IFEM::cout, FieldFuncBase::FIXED_LEVEL, and Real.
Referenced by SIM1D::parse(), SIM2D::parse(), SIM3D::parse(), SIMinput::parse(), SIMinput::parseBCTag(), AnaSol::parseExpressionFunctions(), SIM2D::parseGeometryTag(), SIMoutput::parseOutputTag(), parseRealFunc(), parseTracFunc(), SIMinput::setInitialConditions(), and SIMinput::setNeumann().
| RealFunc * utl::parseRealFunc | ( | const std::string & | func, |
| const std::string & | type = "expression", |
||
| bool | print = true |
||
| ) |
Creates a scalar-valued function by parsing a character string.
| [in] | func | Character string to parse function definition from |
| [in] | type | Function definition type flag |
| [in] | If true, print function definition |
References IFEM::cout, EvalFuncScalar< Scalar >::numError, parseRealFunc(), and Real.
| TensorFunc * utl::parseTensorFunc | ( | const std::string & | func, |
| const std::string & | type | ||
| ) |
Creates a tensor-valued function by parsing a character string.
| [in] | func | Character string to parse function definition from |
| [in] | type | Function definition type flag |
References FieldFuncBase::FIXED_LEVEL, and splitString().
Referenced by AnaSol::parseExpressionFunctions().
| ScalarFunc * utl::parseTimeFunc | ( | const char * | func, |
| const std::string & | type = "expression", |
||
| Real | eps = Real(1.0e-8) |
||
| ) |
Creates a time function by parsing a character string.
| [in] | func | Character string to parse function definition from |
| [in] | type | Function definition type flag |
| [in] | eps | Domain increment for calculation of numerical derivative |
References IFEM::cout.
Referenced by AdaptiveSetup::parse(), SIMCoupledSI< T1, T2 >::parseIterations(), and parseTracFunc().
| TractionFunc * utl::parseTracFunc | ( | const std::string & | func, |
| const std::string & | type = "expression", |
||
| int | dir = 0 |
||
| ) |
Creates a vector-valued function defining a surface traction.
| [in] | func | Character string to parse function definition from |
| [in] | type | Function definition type flag |
| [in] | dir | Coordinate direction of the traction (0=normal direction) |
References IFEM::cout, parseRealFunc(), and Real.
Referenced by SIMinput::parseBCTag(), and SIMinput::setNeumann().
| TractionFunc * utl::parseTracFunc | ( | const tinyxml2::XMLElement * | elem | ) |
Creates a vector-valued function defining a surface traction.
| [in] | elem | Pointer to XML-element to parse function definition from |
References IFEM::cout, getAttribute(), parseRealFunc(), parseTimeFunc(), and parseVecTimeFunc().
| bool utl::parseVec | ( | Vec3 & | value, |
| const char * | argv | ||
| ) |
Parses a character string into a point vector.
| [out] | value | The parsed point vector |
| [in] | argv | Character string with data |
References Real.
| VecFunc * utl::parseVecFunc | ( | const std::string & | func, |
| const std::string & | type = "expression", |
||
| const std::string & | variables = "" |
||
| ) |
Creates a vector-valued function by parsing a character string.
| [in] | func | Character string to parse function definition from |
| [in] | type | Function definition type flag |
| [in] | variables | Optional variable definition for expression functions |
References IFEM::cout, FieldFuncBase::FIXED_LEVEL, EvalFuncScalar< Scalar >::numError, and splitString().
Referenced by AnaSol::parseExpressionFunctions(), and SIMinput::setNeumann().
| VecTimeFunc * utl::parseVecTimeFunc | ( | const char * | func, |
| const std::string & | type | ||
| ) |
Creates a vector-valued time function by parsing a char string.
| [in] | func | Character string to parse function definition from |
| [in] | type | Function definition type flag |
References IFEM::cout.
Referenced by parseTracFunc().
| size_t utl::Pascal | ( | int | p, |
| unsigned short int | nsd | ||
| ) |
Returns the number of monomials in Pascal's triangle.
| [in] | p | Polynomial order (>= 0) |
| [in] | nsd | Number of spatial dimensions (2 or 3) |
Referenced by Pascal().
| char * utl::readLine | ( | std::istream & | is | ) |
Reads one line, ignoring comment lines and leading blanks.
The data read is kept in an internal static buffer.
| is | File stream to read from |
Referenced by AnaSol::AnaSol(), AdaptiveSetup::parse(), NonLinSIM::parse(), SIM1D::parse(), SIM2D::parse(), SIM3D::parse(), SIMinput::parse(), SIMoutput::parse(), TimeStep::parse(), SIMadmin::readFlat(), and ASM::readMatlab().
| bool utl::renumber | ( | int & | num, |
| const IntMap & | old2new, | ||
| bool | msg = false |
||
| ) |
Transforms the integer value num according to a given mapping.
In this method the old2new mapping is not updated.
| num | The integer value to transform, updated value on output | |
| [in] | old2new | Mapping from old values to new values |
| [in] | msg | If true, give error message if given value not found |
| bool utl::renumber | ( | int & | num, |
| int & | runner, | ||
| IntMap & | old2new | ||
| ) |
Transforms the integer value num into a unique range.
This method is invoked on a series of (non-unique) values. When a value that has not been encountered before is detected, it is transformed into runner+1 and runner is incremented.
| num | The integer value to transform, updated value on output |
| runner | The last new value assigned |
| old2new | Mapping from old values to new values return true if the value of num changed |
Referenced by DomainDecomposition::calcGlobalEqNumbers(), DomainDecomposition::calcGlobalNodeNumbers(), GlobalNodes::calcGlobalNodes(), ASMbase::renumberNodes(), and MPC::renumberNodes().
Solves the linear system of equations \( {\bf A x} = {\bf b} \).
If iPivot is null, the equation system is solved directly using DGESV from LAPACK. Otherwise, the matrix A is first LU-factorized using DGETRF (if iPivot[0] is zero), and then the equation system is solved using DGETRS. Therefore, it is possible to invoke this method several times to solve for more right-hand sides, if A and iPivot are retained between each call.
References utl::matrix< T >::cols(), dgesv_(), dgetrf_(), dgetrs_(), M, utl::matrixBase< T >::ptr(), Real, and utl::matrix< T >::rows().
Congruence transformation of a symmetric matrix.
| A | The matrix to be transformed | |
| [in] | Tn | Nodal transformation matrix |
The following matrix multiplication is performed by this function:
\[ {\bf A} = {\bf T}{\bf A}{\bf T}^T \]
where A is a full, symmetric matrix, and the transformation matrix T has the nodal sub-matrix Tn repeated on the diagonal and otherwise zero.
References utl::matrix< T >::cols(), utl::matrix< T >::fill(), K, M, Real, and utl::matrix< T >::rows().
Referenced by SIMsupel::recoverInternalDOFs(), and ASMsupel::transform().
Congruence transformation of a vector.
| V | The vector to be transformed | |
| [in] | Tn | Nodal transformation matrix |
| [in] | transpose | If true, the transpose of Tn is used instead |
The vector V is pre-multiplied with the transformation matrix T which has the nodal sub-matrix Tn (or its transpose) repeated on the diagonal and otherwise zero.
References utl::matrix< T >::cols(), utl::vector< T >::fill(), K, M, Real, utl::matrix< T >::rows(), and utl::vector< T >::size().
|
inline |
Truncate a value to zero when it is less than a given threshold.
Used when printing matrices for easy comparison with other matrices when they contain terms that are numerically zero, except for some round-off noise. The value of the global variable zero_print_tol is used as a tolerance in this method.
References zero_print_tol.
Referenced by SIMoutput::dumpPrimSol(), SIMoutput::dumpResults(), SIMoutput::dumpSolution(), SIMoutput::dumpVector(), NewmarkSIM::initAcc(), operator<<(), SIMoutput::printNRforces(), SIMbase::printSolutionSummary(), SAM::printVector(), NewmarkSIM::solutionNorms(), SIMbase::systemModes(), writeMatlab(), and NonLinSIM::~NonLinSIM().