artisynth.core.femmodels
Class QuadwedgeElement

java.lang.Object
  extended by artisynth.core.modelbase.ModelComponentBase
      extended by artisynth.core.modelbase.RenderableComponentBase
          extended by artisynth.core.femmodels.FemElement
              extended by artisynth.core.femmodels.FemElement3d
                  extended by artisynth.core.femmodels.QuadwedgeElement
All Implemented Interfaces:
CopyableComponent, ModelComponent, RenderableComponent, PropertyChangeListener, ScalableUnits, java.lang.Cloneable, Boundable, HasProperties, HierarchyNode, GLRenderable, GLSelectable, HasRenderProps, Renderable, Scannable

public class QuadwedgeElement
extends FemElement3d


Nested Class Summary
 
Nested classes/interfaces inherited from interface artisynth.core.modelbase.ModelComponent
ModelComponent.NavpanelVisibility
 
Field Summary
 
Fields inherited from class artisynth.core.femmodels.FemElement3d
myProps
 
Fields inherited from class artisynth.core.modelbase.ModelComponentBase
enforceUniqueCompositeNames, enforceUniqueNames, myNumber, NULL_OBJ, useCompactPathNames
 
Fields inherited from interface artisynth.core.modelbase.CopyableComponent
COPY_REFERENCES
 
Fields inherited from interface maspack.render.GLRenderable
TRANSLUCENT, TWO_DIMENSIONAL
 
Constructor Summary
QuadwedgeElement()
           
QuadwedgeElement(FemNode3d[] nodes)
           
QuadwedgeElement(WedgeElement wedge, FemNode3d[] quadraticNodes)
          Create a QuadwedgeElement based on the NODE POSITIONS of a given WedgeElement (i.e.
 
Method Summary
 boolean coordsAreInside(Vector3d coords)
           
 void getdNds(Vector3d dNds, int i, Vector3d coords)
           
 int[] getEdgeIndices()
           
 int[] getFaceIndices()
           
 double getH(int i, Vector3d coords)
          Returns the value of the pressure shape function.
 double[] getIntegrationCoords()
           
 IntegrationPoint3d[] getIntegrationPoints()
           
 double getN(int i, Vector3d coords)
           
 double[] getNodalExtrapolationMatrix()
           
 double[] getNodeCoords()
           
 Matrix getPressureWeightMatrix()
          Returns the pressure weight matrix for this element.
static FemNode3d[] getQuadraticNodes(FemNode3d n0, FemNode3d n1, FemNode3d n2, FemNode3d n3, FemNode3d n4, FemNode3d n5)
           
static FemNode3d[] getQuadraticNodes(WedgeElement wedge)
           
 IntegrationPoint3d getWarpingPoint()
           
 boolean isInside(Point3d pnt)
          Tests whether or not a point is inside an element.
 int numIntegrationPoints()
           
 int numPressureVals()
          Returns the number of pressure variables associated with this element.
 void renderWidget(GLRenderer renderer, double size, RenderProps props)
           
 FemNode3d[][] triangulateFace(FaceNodes3d face)
           
 
Methods inherited from class artisynth.core.femmodels.FemElement3d
addAuxiliaryMaterial, addNodeForce, addNodeForce0, addNodeForce0, addNodeStiffness, addNodeStiffness, clearState, computeCentroid, computeCovariance, computeDirectedRenderSize, computeGravityWeights, computePressures, computeRenderCoordsAndGradient, computeRestVolumes, computeVolumes, computeWarping, computeWarping, connectToHierarchy, copy, createElement, createElement, createIntegrationPoints, disconnectFromHierarchy, getAllPropertyInfo, getAuxiliaryMaterials, getCopyReferences, getElementWidgetSize, getElementWidgetSizeMode, getFrame, getIncompressConstraints, getIncompressIndex, getIntegrationData, getMarkerCoordinates, getNaturalCoordinates, getNaturalCoordinatesRobust, getNodeCoords, getNodes, getNumEdges, getNumFaces, getPoint, getTriFaces, getWarpingData, hasEdge, hasFace, hasFace, invalidateRestData, isInvertedAtRest, materialsAreInvertible, numAuxiliaryMaterials, numPoints, removeAuxiliaryMaterial, render, render, renderEdges, renderWidget, renderWidget, renderWidgetEdges, renderWidgetEdges, scaleDistance, setElementWidgetSize, setElementWidgetSizeMode, setFrame, setIncompressIndex, updateBounds, updateWarpingStiffness
 
Methods inherited from class artisynth.core.femmodels.FemElement
containsNode, createRenderProps, getDensity, getDensityMode, getEffectiveMaterial, getHardReferences, getIndex, getLocalNodeIndex, getMass, getMaterial, getRestVolume, getSelection, getVolume, hasActiveNodes, hasControllableNodes, integrationPointsMapToNodes, isDuplicatable, isInverted, numberString, numNodes, prerender, propertyChanged, scaleMass, setDensity, setDensityMode, setIndex, setInverted, setMass, setMaterial, updateNodeMasses
 
Methods inherited from class artisynth.core.modelbase.RenderableComponentBase
getRenderHints, getRenderProps, isSelectable, numSelectionQueriesNeeded, setRenderProps, updateRenderProps
 
Methods inherited from class artisynth.core.modelbase.ModelComponentBase
checkFlag, checkName, checkNameUniqueness, clearFlag, clone, createTempFlag, getChildren, getGrandParent, getName, getNameRange, getNavpanelVisibility, getNavpanelVisibility, getNumber, getParent, getProperty, getSoftReferences, hasChildren, hasState, isFixed, isMarked, isSelected, isWritable, makeValidName, makeValidName, notifyParentOfChange, postscan, printReferences, recursivelyContained, recursivelyContains, removeTempFlag, scan, setFixed, setFlag, setMarked, setName, setNavpanelVisibility, setNavpanelVisibility, setNumber, setParent, setSelected, updateReferences, write
 
Methods inherited from class java.lang.Object
equals, getClass, hashCode, notify, notifyAll, toString, wait, wait, wait
 
Methods inherited from interface artisynth.core.modelbase.ModelComponent
getName, getNavpanelVisibility, getNumber, getParent, getSoftReferences, hasState, isFixed, isMarked, isSelected, notifyParentOfChange, postscan, scan, setFixed, setMarked, setName, setNumber, setParent, setSelected, updateReferences
 
Methods inherited from interface maspack.properties.HasProperties
getProperty
 
Methods inherited from interface maspack.properties.HierarchyNode
getChildren, hasChildren
 
Methods inherited from interface maspack.util.Scannable
isWritable, write
 

Constructor Detail

QuadwedgeElement

public QuadwedgeElement()

QuadwedgeElement

public QuadwedgeElement(WedgeElement wedge,
                        FemNode3d[] quadraticNodes)
Create a QuadwedgeElement based on the NODE POSITIONS of a given WedgeElement (i.e. it does not inherit any other attributes of the WedgeElement). Takes the 6 nodes of the given WedgeElement along with the 9 given quadraticNodes as the middle nodes to create a 15 node quadratic wedge.

Parameters:
wedge - A wedge element
quadraticNodes - the 9 edge nodes

QuadwedgeElement

public QuadwedgeElement(FemNode3d[] nodes)
Method Detail

getIntegrationPoints

public IntegrationPoint3d[] getIntegrationPoints()
Specified by:
getIntegrationPoints in class FemElement3d

getWarpingPoint

public IntegrationPoint3d getWarpingPoint()
Specified by:
getWarpingPoint in class FemElement3d

coordsAreInside

public boolean coordsAreInside(Vector3d coords)
Specified by:
coordsAreInside in class FemElement3d

numIntegrationPoints

public int numIntegrationPoints()
Specified by:
numIntegrationPoints in class FemElement3d

getIntegrationCoords

public double[] getIntegrationCoords()
Specified by:
getIntegrationCoords in class FemElement3d

getNodeCoords

public double[] getNodeCoords()
Specified by:
getNodeCoords in class FemElement3d

getNodalExtrapolationMatrix

public double[] getNodalExtrapolationMatrix()
Specified by:
getNodalExtrapolationMatrix in class FemElement3d

getN

public double getN(int i,
                   Vector3d coords)
Specified by:
getN in class FemElement3d

numPressureVals

public int numPressureVals()
Returns the number of pressure variables associated with this element. All of the linear elements have one pressure variable, whereas some of the quadratic elements have more.

Overrides:
numPressureVals in class FemElement3d
Returns:
number of pressure variables.make

getH

public double getH(int i,
                   Vector3d coords)
Returns the value of the pressure shape function. By default, this method returns 1, corresponding to a single pressure variable with constant value over the entire element. Elements with a larger number of pressure DOFs should override this method to supply the appropriate shape functions.

Overrides:
getH in class FemElement3d
Parameters:
i - index of the pressure variable; should be less than the value returned by FemElement3d.numPressureVals()
coords - coordinates at which the shape function should be evaluated.

getPressureWeightMatrix

public Matrix getPressureWeightMatrix()
Description copied from class: FemElement3d
Returns the pressure weight matrix for this element. The pressure weight matrix is given by the inverse of the integral of H^T H, where H is the row vector formed from the pressure shape functions.

By default, this method returns a pressure weight matrix for the case where there is only one pressure value. Such matrices always have a single value of 1. Elements with a larger number of pressure values should override this method to return a pressure weight matrix appropriate for that element.

Overrides:
getPressureWeightMatrix in class FemElement3d

getdNds

public void getdNds(Vector3d dNds,
                    int i,
                    Vector3d coords)
Specified by:
getdNds in class FemElement3d

getQuadraticNodes

public static FemNode3d[] getQuadraticNodes(WedgeElement wedge)

getQuadraticNodes

public static FemNode3d[] getQuadraticNodes(FemNode3d n0,
                                            FemNode3d n1,
                                            FemNode3d n2,
                                            FemNode3d n3,
                                            FemNode3d n4,
                                            FemNode3d n5)

getEdgeIndices

public int[] getEdgeIndices()
Specified by:
getEdgeIndices in class FemElement3d

getFaceIndices

public int[] getFaceIndices()
Specified by:
getFaceIndices in class FemElement3d

renderWidget

public void renderWidget(GLRenderer renderer,
                         double size,
                         RenderProps props)
Overrides:
renderWidget in class FemElement3d

triangulateFace

public FemNode3d[][] triangulateFace(FaceNodes3d face)
Overrides:
triangulateFace in class FemElement3d

isInside

public boolean isInside(Point3d pnt)
Description copied from class: FemElement3d
Tests whether or not a point is inside an element.

Overrides:
isInside in class FemElement3d
Parameters:
pnt - point to check if is inside
Returns:
true if point is inside the element