public class MechSystemSolver
extends java.lang.Object
Modifier and Type | Class and Description |
---|---|
static class |
MechSystemSolver.Integrator |
static class |
MechSystemSolver.PosStabilization
Indicates the method by which positions should be stabilized.
|
Modifier and Type | Field and Description |
---|---|
boolean |
alwaysProjectFriction |
static boolean |
DEFAULT_HYBRID_SOLVES_ENABLED |
static boolean |
myAlwaysAnalyze |
static boolean |
myDefaultHybridSolveP |
boolean |
myMurtyVelSolveRebuild |
static double |
myT1 |
boolean |
myUseImplicitFriction |
static int |
NO_SYS_UPDATE
Flag for KKTSolverFactorAndSolve methods: do not update the system with
the computed velocities and constraint forces.
|
boolean |
profileConstrainedBE |
static boolean |
profileConstraintSolves |
boolean |
profileImplicitFriction |
boolean |
profileKKTSolveTime |
boolean |
profileWholeSolve |
static int |
TRAPEZOIDAL
Flag for KKTSolverFactorAndSolve methods: indicates trapezoidal
integration
|
static boolean |
useFictitousJacobianForces |
Constructor and Description |
---|
MechSystemSolver(MechSystem system)
Create a new MechSystem solver for a specified MechSystem.
|
MechSystemSolver(MechSystem system,
MechSystemSolver solver)
Create a new MechSystem solver for a specified MechSystem,
with settings imported from an existing MechSystemSolver.
|
Modifier and Type | Method and Description |
---|---|
void |
addActiveMassMatrix(MechSystem sys,
SparseBlockMatrix S) |
void |
addMassForces(VectorNd f,
double t) |
void |
addSubVector(VectorNd vr,
VectorNd v1,
int off,
int size)
Computes
|
void |
backwardEuler(double t0,
double t1,
StepAdjustment stepAdjust) |
static double |
brentRootFinder(double a,
double fa,
double b,
double fb,
double eps,
double feps,
Function1x1 func)
Implementation of Brent's root-finding method, guarantees that the final call
to func is at the root
|
void |
computeParametricForces(double h) |
void |
constrainedBackwardEuler(double t0,
double t1,
StepAdjustment stepAdjust) |
void |
constrainedForwardEuler(double t0,
double t1,
StepAdjustment stepAdjust) |
SparseBlockMatrix |
createActiveBilateralMatrix(double t) |
SparseBlockMatrix |
createActiveDampingMatrix(double h) |
SparseBlockMatrix |
createActiveMassMatrix(double t)
Creates and returns a sparse block matrix consisting of the current
active portion of the mass matrix.
|
SparseBlockMatrix |
createActiveStiffnessMatrix(double h) |
void |
dispose() |
void |
finalize() |
void |
fullBackwardEuler(double t0,
double t1,
StepAdjustment stepAdjust) |
static boolean |
getAlwaysAnalyze() |
int |
getFrictionIterations() |
boolean |
getHybridSolve() |
static boolean |
getHybridSolvesEnabled()
Queries whether hybrid sparse matrix solves are enabled by default.
|
MechSystemSolver.Integrator |
getIntegrator() |
VectorNd |
getLambda() |
SparseSolverId |
getMatrixSolver() |
int |
getMaxIterations() |
SparseBlockMatrix |
getSolveMatrix()
Returns the solve matrix.
|
MechSystemSolver.PosStabilization |
getStabilization() |
int |
getStaticIncrements()
Returns the number of load increments being used with the
MechSystemSolver.Integrator.StaticIncremental
integrator. |
double |
getStaticTikhonovFactor() |
double |
getTolerance() |
IterativeSolver.ToleranceType |
getToleranceType() |
boolean |
getUpdateForcesAtStepEnd() |
boolean |
getUseImplicitFriction() |
void |
initialize() |
void |
KKTFactorAndSolve(VectorNd vel,
VectorNd fpar,
VectorNd bf,
VectorNd btmp,
VectorNd vel0,
double h,
double a0,
double a1,
double a2,
double a3,
int flags)
Solves a KKT system in which the Jacobian augmented M matrix and
and force vectors are given by
|
void |
KKTFactorAndSolve(VectorNd vel,
VectorNd fpar,
VectorNd bf,
VectorNd btmp,
VectorNd vel0,
double h,
int flags)
Solves a KKT system in which the Jacobian augmented M matrix and
and force vectors are given by
|
void |
KKTSolve(VectorNd vel,
VectorNd lam,
VectorNd the,
VectorNd bf)
Solves the KKT system given a new right hand side
bf . |
void |
KKTStaticFactorAndSolve(VectorNd u,
VectorNd bf,
double beta,
VectorNd btmp)
Solves a static KKT system of the form
|
static double |
modifiedGoldenSection(double a,
double fa,
double b,
double fb,
double eps,
double feps,
Function1x1 func)
Implementation of a modified Golden section search for minimizing |f(s)|,
guaranteeing that the final call to func is at the returned minimum
|
void |
mulActiveInertias(VectorNd b,
VectorNd v) |
void |
nonDynamicSolve(double t0,
double t1,
StepAdjustment stepAdjust) |
void |
projectPosConstraints(double t) |
void |
projectRigidBodyPosConstraints(double t) |
static void |
setAlwaysAnalyze(boolean enable) |
void |
setCrsFileName(java.lang.String name) |
void |
setFrictionIterations(int num) |
void |
setHybridSolve(boolean enable) |
static void |
setHybridSolvesEnabled(boolean enable)
Set hybrid sparse matrix solves to be enabled or disabled by default.
|
void |
setIntegrator(MechSystemSolver.Integrator integrator) |
void |
setIterativeSolver(IterativeSolver solver) |
static void |
setLogWriter(java.io.PrintWriter writer) |
void |
setMatrixSolver(SparseSolverId solver) |
void |
setMaxIterations(int max) |
void |
setParametricTargets(double s,
double h) |
void |
setStabilization(MechSystemSolver.PosStabilization stabilization) |
void |
setStaticIncrements(int iters)
Sets the number of load increments to use with the
MechSystemSolver.Integrator.StaticIncremental
integrator. |
void |
setStaticTikhonovFactor(double eps)
Applies Tiknonov regularization for static solves, minimizing
|
void |
setSubVector(VectorNd vr,
VectorNd v1,
int off,
int size)
Sets
|
void |
setTolerance(double tol) |
void |
setToleranceType(IterativeSolver.ToleranceType type) |
void |
setUpdateForcesAtStepEnd(boolean enable) |
void |
setUseImplicitFriction(boolean enable) |
void |
solve(double t0,
double t1,
StepAdjustment stepAdjust) |
void |
staticIncremental(double t0,
double t1,
int nincrements,
StepAdjustment stepAdjust) |
void |
staticIncrementalStep(double t0,
double t1,
double alpha,
StepAdjustment stepAdjust)
Scales forces and constraints down by alpha, and solves the adjusted problem
|
void |
staticLineSearch(double t0,
double t1,
StepAdjustment stepAdjust) |
void |
trapezoidal(double t0,
double t1,
StepAdjustment stepAdjust) |
void |
updateConstraintMatrices(double h,
boolean includeFriction) |
void |
updateMassMatrix(double t) |
void |
updateStateSizes() |
boolean |
usingImplicitFriction() |
public static int NO_SYS_UPDATE
public static int TRAPEZOIDAL
public boolean profileKKTSolveTime
public boolean profileWholeSolve
public boolean profileConstrainedBE
public boolean profileImplicitFriction
public boolean alwaysProjectFriction
public boolean myMurtyVelSolveRebuild
public boolean myUseImplicitFriction
public static boolean useFictitousJacobianForces
public static boolean myAlwaysAnalyze
public static boolean DEFAULT_HYBRID_SOLVES_ENABLED
public static boolean myDefaultHybridSolveP
public static boolean profileConstraintSolves
public static double myT1
public MechSystemSolver(MechSystem system)
public MechSystemSolver(MechSystem system, MechSystemSolver solver)
public void setUpdateForcesAtStepEnd(boolean enable)
public boolean getUpdateForcesAtStepEnd()
public static boolean getAlwaysAnalyze()
public static void setAlwaysAnalyze(boolean enable)
public static boolean getHybridSolvesEnabled()
true
if hybrid solves are enabledpublic static void setHybridSolvesEnabled(boolean enable)
enable
- if true
, enables hybrid solvespublic boolean getHybridSolve()
public void setHybridSolve(boolean enable)
public boolean getUseImplicitFriction()
public void setUseImplicitFriction(boolean enable)
public boolean usingImplicitFriction()
public int getFrictionIterations()
public void setFrictionIterations(int num)
public void setParametricTargets(double s, double h)
public void computeParametricForces(double h)
public void updateMassMatrix(double t)
public SparseBlockMatrix createActiveMassMatrix(double t)
t
- current timepublic void updateStateSizes()
public void setIntegrator(MechSystemSolver.Integrator integrator)
public void setIterativeSolver(IterativeSolver solver)
public MechSystemSolver.Integrator getIntegrator()
public double getTolerance()
public void setTolerance(double tol)
public IterativeSolver.ToleranceType getToleranceType()
public void setToleranceType(IterativeSolver.ToleranceType type)
public void setMaxIterations(int max)
public int getMaxIterations()
public void setMatrixSolver(SparseSolverId solver)
public void setStabilization(MechSystemSolver.PosStabilization stabilization)
public MechSystemSolver.PosStabilization getStabilization()
public SparseSolverId getMatrixSolver()
public void nonDynamicSolve(double t0, double t1, StepAdjustment stepAdjust)
public void solve(double t0, double t1, StepAdjustment stepAdjust)
public VectorNd getLambda()
public void addActiveMassMatrix(MechSystem sys, SparseBlockMatrix S)
public void addMassForces(VectorNd f, double t)
public void addSubVector(VectorNd vr, VectorNd v1, int off, int size)
vr += v1(off:off+size-1)
public void setSubVector(VectorNd vr, VectorNd v1, int off, int size)
vr = v1(off:off+size-1)
public void backwardEuler(double t0, double t1, StepAdjustment stepAdjust)
public void updateConstraintMatrices(double h, boolean includeFriction)
public void KKTFactorAndSolve(VectorNd vel, VectorNd fpar, VectorNd bf, VectorNd btmp, VectorNd vel0, double h, int flags)
M' = M - h df/dv - h^2 df/dx bf' = bf + (-h df/dv) vel0When used to solve for velocities in an implicit integrator, then on input, bf is assumed to be given by
bf = M vel0 + h fwhere h is the time step and f is the generalized forces, while on output bf is modified to include the Jacobian terms described above. XXX question about which time do we want? t0 or t1? Right now it is set to t1; setting it to t0 seems to cause a small difference in the inverse tongue model.
vel
- returns the computed velocityfpar
- if useFictitousJacobianForces is true, returns fictitious
Jacobian forces for parametric componentsbf
- right side offsetbtmp
- temporary vectorvel0
- right side velocityh
- interval time stepflags
- NO_SYS_UPDATEpublic void KKTFactorAndSolve(VectorNd vel, VectorNd fpar, VectorNd bf, VectorNd btmp, VectorNd vel0, double h, double a0, double a1, double a2, double a3, int flags)
M' = M + a0 df/dv + a1 df/dx bf' = bf + (a2 df/dv + a3 df/dx) vel0It is assumed that a0 and a1 are both non-zero. It is also assumed that the a0 = -alpha h, where h is the step size and alpha indicates the propertion of implicitness for the solve; i.e., for regular backward euler, alpha=1, while for trapezoidal solves, alpha = 0.5; When used to solve for velocities in an implicit integrator, then on input, bf is assumed to be given by
bf = M vel0 + h fwhere h is the time step and f is the generalized forces, while on output bf is modified to include the Jacobian terms described above.
vel
- returns the computed velocityfpar
- if useFictitousJacobianForces is true, returns fictitious
Jacobian forces for parametric componentsbf
- right side offsetbtmp
- temporary vectorvel0
- right side velocityh
- interval time step - used to scale constraint offsets and
impulsesa0
- left side df/dv coefficienta1
- left side df/dx coefficienta2
- right side df/dv coefficienta3
- right side df/dx coefficientflags
- NO_SYS_UPDATE, TRAPEZOIDALpublic void KKTStaticFactorAndSolve(VectorNd u, VectorNd bf, double beta, VectorNd btmp)
-df/dx*Delta(x) -G^T*lambda - N^T*theta = f
G*Delta(x) + g = 0, N*Delta(x) + n >= 0
u
- returned displacement Delta(x)bf
- right-hand side net forcebeta
- scale factor for any additional forces such as fictitious forcesbtmp
- temporary vectorpublic void setCrsFileName(java.lang.String name)
public static void setLogWriter(java.io.PrintWriter writer)
public void KKTSolve(VectorNd vel, VectorNd lam, VectorNd the, VectorNd bf)
bf
. It is
assumed that KKTFactorAndSolve() has already been called once, and that
the initial velocity vel0
is unchanged from that
call. The updated value for bf
should be based
on the value returned by KKTFactorAndSolve(), as described below.
Recall that on input to KKTFactorAndSolve(), bf
is
assumed to be of the form
bf = M vel0 + h f1where
M
is the mass matrix, h
is the time step,
and f1
are the force values. On return, bf
is
modified by adding values that depend on vel0
and the
position and velocity Jacobians of the system. These additional values
are needed on subsequent calls to KKTSolve(). Therefore, when using
KKTSolve() to compute velocities for different force values
f2
, one should use the bf
returned by
KKTFactorAndSolve() and simply adjust it to reflect f2
instead of f1
:
bf1 = M vel0 + h f1 KKTFactorAndSolve (vel1, bf1, btmp, vel0, ...) bf2 = bf1 - h f1 + h f2 KKTSolve (vel2, myLam, the, bf2)
public void projectPosConstraints(double t)
public void projectRigidBodyPosConstraints(double t)
public void constrainedForwardEuler(double t0, double t1, StepAdjustment stepAdjust)
public void constrainedBackwardEuler(double t0, double t1, StepAdjustment stepAdjust)
public void fullBackwardEuler(double t0, double t1, StepAdjustment stepAdjust)
public void trapezoidal(double t0, double t1, StepAdjustment stepAdjust)
public void setStaticTikhonovFactor(double eps)
L(dx) = W(x+dx) + eps*|dx|^2where W is the energy potential function (including constraints) for the static system. If eps is less than zero, then we try to estimate the optimal parameter based on the Frobenius norm of the stiffness matrix.
eps
- tikhonov regularization factorpublic double getStaticTikhonovFactor()
public void setStaticIncrements(int iters)
MechSystemSolver.Integrator.StaticIncremental
integrator.iters
- number of load incrementspublic int getStaticIncrements()
MechSystemSolver.Integrator.StaticIncremental
integrator.MechSystemSolver.Integrator.StaticIncremental
public void staticIncrementalStep(double t0, double t1, double alpha, StepAdjustment stepAdjust)
t0
- start time for solvet1
- stop time for solvealpha
- step factorstepAdjust
- step adjustment descriptionpublic void staticIncremental(double t0, double t1, int nincrements, StepAdjustment stepAdjust)
public static double brentRootFinder(double a, double fa, double b, double fb, double eps, double feps, Function1x1 func)
a
- left-side of intervalfa
- function evaluated at ab
- right-side of intervalfb
- function evaluated at beps
- tolerance for interval [a,b]feps
- tolerance for function value |f| <
feps considered rootfunc
- function to evaluatepublic static double modifiedGoldenSection(double a, double fa, double b, double fb, double eps, double feps, Function1x1 func)
a
- left-side of search intervalfa
- f(a)b
- right-side of search intervalfb
- f(b)eps
- tolerance for interval [a,b]feps
- tolerance for function evaluation, |f(s)| <
feps considered a rootfunc
- function to evaluatepublic void staticLineSearch(double t0, double t1, StepAdjustment stepAdjust)
public SparseBlockMatrix createActiveStiffnessMatrix(double h)
public SparseBlockMatrix createActiveDampingMatrix(double h)
public SparseBlockMatrix createActiveBilateralMatrix(double t)
public SparseBlockMatrix getSolveMatrix()
public void dispose()
public void finalize()
finalize
in class java.lang.Object
public void initialize()