gaussMatch: map two surface meshes using smoothed displacement fields

View source: R/gaussDisplace.r

gaussMatchR Documentation

map two surface meshes using smoothed displacement fields

Description

Map a reference mesh onto a target surface using displacement fields.

Usage

gaussMatch(
  x,
  mesh2,
  iterations = 10,
  smooth = NULL,
  smoothit = 10,
  smoothtype = c("taubin", "laplace", "HClaplace"),
  sigma = 20,
  displacementsmooth = c("Gauss", "Laplace", "Exponential"),
  gamma = 2,
  f = 1.2,
  oneway = F,
  lm1 = NULL,
  lm2 = NULL,
  rigid = NULL,
  similarity = NULL,
  affine = NULL,
  tps = FALSE,
  nh = NULL,
  toldist = 0,
  pro = c("kd", "vcg", "morpho"),
  k0 = 50,
  prometh = 1,
  angtol = pi/2,
  border = FALSE,
  horiz.disp = NULL,
  useiter = FALSE,
  AmbergK = NULL,
  AmbergLambda = NULL,
  tol = 1e-05,
  useConstrained = TRUE,
  angclost = TRUE,
  noinc = FALSE,
  silent = FALSE,
  visualize = FALSE,
  folder = NULL,
  alpha = 0.7,
  col1 = "red",
  col2 = "white",
  add = FALSE,
  bbox = NULL,
  bboxCrop = NULL,
  threads = 0,
  cb = NULL,
  useValid2Constrain = FALSE,
  mahasafe = 1e+10,
  ...
)

Arguments

x

reference mesh: triangular mesh of class "mesh3d"or of class BayesDeform created by createBayes to restrict based on a known distribution. To use this option the package RvtkStatismo https://github.com/zarquon42b/RvtkStatismo has to be installed. If x is a model, it works best if mesh2 is already aligned to the model's mean.

mesh2

An object of class mesh3d used as target mesh. Mesh resolution should be ~1.5.

iterations

Iterations of displacement. Default is 10.

smooth

Integer: smooth the resulting mesh at the smooth-th iteration. E.g. a value of 2 will smooth the resulting mesh at each 2nd iteration. Not to be confused with displacementsmooth!! Default is NULL, no smoothing.

smoothit

integer: smoothing steps.

smoothtype

Type of smoothing: Taubin, Laplacian, or HClaplacian. For details see vcgSmooth.

sigma

starting value of kernel bandwidth/B-spline support. For all kernels except B-spline, sigma controls the importance of the neighbourhood by defining the bandwidth of the smoothing kernel. For B-spline it defines the support (the higher, the "wobblier" the deformation field can become.

displacementsmooth

kernel function for smoothing are "Gauss","Laplace", "Exponential" and "Bspline" (or any abbreviation thereof).

gamma

dampening factor controlling displacement strength. The smoothed displacement vector for each vertex is divided by gamma. The larger gamma, the slower the approximation.

f

parameter controlling iterative decrease (if > 1, increase else) of sigma making the displacement locally more elastic with each iteration. Starting with sigma, this parameter for the k-th iteration is sigma *f ^(-k)

oneway

logical: only displace towards the target without taking into account the displacement from the target.

lm1

A k x 3 matrix containing landmarks corrresponding to mesh1 for initial rotation of mesh1 onto mesh2.

lm2

A k x 3 matrix containing landmarks corrresponding to mesh2 for initial rotation of mesh1 onto mesh2.

rigid

named list. Passing parameters to icp, for rigid registration. If landmarks are provided and only those should count, set rigid$iterations=0.

similarity

named list. Passing parameters to icp, for similarity registration (rigid +scaling). If landmarks are provided and only those should count, set similarity$iterations=0 (and rigid=NULL).

affine

named list. Passing parameters to icp, for affine registration. If landmarks are provided and only those should count, set similarity$iterations=0 (with rigid=NULL and similarity=NULL)

tps

logical: if TRUE and landmarks are provided, the reference will be mapped to the target using a Thin-Plate Spline interpolation. Overrides rigid,affine and similarity.

nh

Integer: neighbourhood (number vertices) for controlling displacement smoothing, default is 150/mesh resolution.

toldist

Integer: Exclude everything from the whole procedure with a greater distance from initial point than toldist. 0 disables this feature.

pro

which projection method to use: "kd" = vcgClostKD, "m"= closemeshKD from Morpho; "v"= vcgClost from package Rvcg.

k0

Integer: argument passed to closemeshKD (will be argument "k" in closemeshKD .

prometh

argument passed to closemeshKD. Integer: 0 or 1. If prometh=0, take closest point for displacement. If prometh=1, do not just take the closest point, but for two absolut distances which are the same, take the point which is orthogonal to the closest face see Moshfeghi 1994).

angtol

numeric: If the angle between hit points' normals and the starting points' normals exceeds this threshold the displacement vector will be discarded.

border

Logical: if TRUE, displacement vectors hitting mesh borders are discarded.

horiz.disp

numeric: If the angle between hit points' normals (independent of its orientation) and the distance vector between hit point and starting points exceeds this threshold, the displacement vector will be discarded. Reduces distortion especially at mesh borders.

useiter

logical: if AmbergK and AmbergLambda are set and useiter is TRUE, the minimization will be performed using the latest iteration (slower). The orginal reference otherwise.

AmbergK

a single integer or an integer vector vector containing the k0-value (normal slackness) for each iteration for a smooth Deformation using AmbergDeformSpam.

AmbergLambda

as single numeric value or a numeric vector containing the lambda-value for each iteration for a smooth Deformation using AmbergDeformSpam.

tol

convergence threshold: if RMSE between iterations is below tol, the function stops.

useConstrained

logical: if TRUE and Bayes and landmarks are defined, the landmarks are not only used to get a suitable reference but the model will also be constrained by the landmarks to subsequently restrict the shape variability. If FALSE, the full model is used.

angclost

if TRUE, the closest k faces will be evaluated and the closest with the appropriate normal angle will be selected.

noinc

after each iteration the RMSE between target and moving image is calculated and if this value increases compared to a previous value, the matching stops. Can be useful when matching a statistical model to a partial shape.

silent

logical suppress messages

visualize

logical: if TRUE the matching process is visualized

folder

character: if specified, each a screenshot of each deformation step will be saved as a png file in this folder.

alpha

numeric between 0 and 1 controls opacity of target mesh if visualize=TRUE.

col1

color of fix mesh (if visualize = TRUE)

col2

color of moving mesh (if visualize = TRUE)

add

logical: if FALSE, the 3D window will be reset.

bbox

extend of the margins around the target shape to be considered.

bboxCrop

extend of the bounding box around mesh1 (after alignmend) that will be cropped from target to speed things up.

threads

integer: threads to use in multithreaded routines.

cb

optional: callback function that takes arguments i="current iteration", distance= "distance from target to current estimate" and t.dist="average vertex displacement to last iteration"

useValid2Constrain

logical: if TRUE and x is a shape model, then only those vertices with valid hits are used to compute the PosteriorMean.

mahasafe

numeric: define the max allowed per-vertex mahalanobis distance (only available if the fitting is model based). Currently only working with the mahasafe branch of RvtkStatismo.

...

Further arguments passed to nn2.

Details

This function implements the mesh matching method suggested by Moshfeghi et al. and Bryan et al.. Additional mechanisms for controlling and restricting the displacement smoothing are implemented

Value

If a patch is specified:

mesh

matched mesh

patch

displaced patch as specified in input.

else a mesh of class "mesh3d" is returned.

Note

Based on the closest points (constrained by the various options), the displacement field will be smoothed using kernel functions of the k-closest displacement vectors. The smoothing kernels are "Gauss","Laplace", "Exponential" and "Bspline". The displacement at point x will be the weighted displacment vectors of the k-closest displacement vectors. Be d the distance to a neightbouring point, the weight will be calculated as:

Gaussian: w(d) = exp(\frac{-d^2}{2\sigma^2})

Laplacian: w(d) = exp(\frac{-d}{\sigma})

Exponential: w(d) = exp(\frac{-d}{2\sigma^2})

Author(s)

Stefan Schlager

References

Bryan, R., Mohan, P. S., Hopkins, A., Galloway, F., Taylor, M., and Nair, P. B. 2010. Statistical modelling of the whole human femur incorporating geometric and material properties. Medical Engineering & Physics, 32(1):57-65.

Moshfeghi, M., Ranganath, S., and Nawyn, K. 1994. Three-dimensional elastic matching of volumes. IEEE Transactions on Image Processing: A Publication of the IEEE Signal Processing Society, 3(2):128-138.

See Also

meshres, vcgClost, vcgBorder, icp, vcgSmooth

outsideBBox,getMeshBox

Examples

require(Morpho)
data(nose)##load data
##warp a mesh onto another landmark configuration:
longnose.mesh <- tps3d(shortnose.mesh,shortnose.lm,longnose.lm)
### result won't be too good as the surfaces do stronly differ.
## we start with an affine transformation initiated by landmarks
affine <- list(iterations=20,subsample=100,rhotol=pi/2,uprange=0.9)
match <- gaussMatch(shortnose.mesh,longnose.mesh,lm1=shortnose.lm,
                   lm2=longnose.lm,gamma=2,iterations=10,smooth=1,smoothtype="h",
                   smoothit=10,nh=50,angtol=pi/2,affine=affine,sigma=20)

## Not run: 
## now a more cautiously approach using a larger neighbourhood
## (400 instead of 50) and no intermediary smoothing:
matchNoSmooth <- gaussMatch(shortnose.mesh,longnose.mesh,lm1=shortnose.lm,
                   lm2=longnose.lm,gamma=2,iterations=20,,nh=400,angtol=pi/2,
                   affine=affine,sigma=20)

## End(Not run)

zarquon42b/mesheR documentation built on Jan. 28, 2024, 2:17 p.m.