.packageName <- "cluster"
## Ensure consistent "diss.." class --- make "namespace-private-global !
dissiCl <- c("dissimilarity", "dist")

## consistent error / warning messages; could use for internationalization
..msg <-
    list(error =
         c(NAdiss = "NA-values in the dissimilarity matrix not allowed.",
           non.diss="x is not and cannot be converted to class dissimilarity"
           ),
         warn = c()
         )
#### $Id: agnes.q,v 1.17 2004/03/11 16:26:40 maechler Exp maechler $
agnes <- function(x, diss = inherits(x, "dist"), metric = "euclidean",
		  stand = FALSE, method = "average", par.method,
                  keep.diss = n < 100, keep.data = !diss)
{
    METHODS <- c("average", "single","complete", "ward","weighted", "flexible")
    ## hclust has more;  1    2         3           4       5         6
    meth <- pmatch(method, METHODS)
    if(is.na(meth)) stop("invalid clustering method")
    if(meth == -1) stop("ambiguous clustering method")
    method <- METHODS[meth]
    if(method == "flexible") {
        ## Lance-Williams formula (but *constant* coefficients):
        par.method <- as.numeric(par.method) # or barf
        stopifnot((np <- length(par.method)) >= 1)
        if(np == 1)## default (a1= a, a2= a, b= 1-2a, c = 0)
            par.method <- c(par.method, par.method, 1-2*par.method, 0)
        else if(np == 3)
            par.method <- c(par.method, 0)
        else if(np != 4)
            stop("'par.method' must be of length 1, 3, or 4")
        attr(method,"par") <- par.method
    } else par.method <- double(1)

    if((diss <- as.logical(diss))) {
	## check type of input vector
	if(any(is.na(x))) stop(..msg$error["NAdiss"])
	if(data.class(x) != "dissimilarity") { # try to convert to
	    if(!is.null(dim(x))) {
		x <- as.dist(x) # or give an error
	    } else {
		## possibly convert input *vector*
		if(!is.numeric(x) || is.na(n <- sizeDiss(x)))
		    stop(..msg$error["non.diss"])
		attr(x, "Size") <- n
	    }
	    class(x) <- dissiCl
	    if(is.null(attr(x,"Metric"))) attr(x, "Metric") <- "unspecified"
	}
	n <- attr(x, "Size")
	dv <- x[lower.to.upper.tri.inds(n)]
	## prepare arguments for the Fortran call
	dv <- c(0., dv)# "double", 1st elem. "only for Fortran" (?)
	jp <- 1
	mdata <- FALSE
	ndyst <- 0
	x2 <- double(1)
	jdyss <- 1 # distances on _input_
    }
    else {
	## check input matrix and standardize, if necessary
	x <- data.matrix(x)
	if(!is.numeric(x)) stop("x is not a numeric dataframe or matrix.")
	x2 <- if(stand) scale(x, scale = apply(x, 2, meanabsdev)) else x
        storage.mode(x2) <- "double"
	ndyst <- if(metric == "manhattan") 2 else 1
	n <- nrow(x2)
	jp <- ncol(x2)
	if((mdata <- any(inax <- is.na(x2)))) { # TRUE if x[] has any NAs
	    jtmd <- as.integer(ifelse(apply(inax, 2, any), -1, 1))
	    ## VALue for MISsing DATa
	    valmisdat <- 1.1* max(abs(range(x2, na.rm=TRUE)))
	    x2[inax] <- valmisdat
	    valmd <- rep(valmisdat, jp)
	}
	dv <- double(1 + (n * (n - 1))/2)
	jdyss <- 0 # distances to be computed
    }
    if(keep.diss) jdyss <- jdyss + 10
    ## call Fortran routine
    res <- .Fortran("twins",
		    as.integer(n),
		    as.integer(jp),
		    x2,
		    dv,
		    dis = double(if(keep.diss) length(dv) else 1),
		    ok = as.integer(jdyss),
		    if(mdata) valmd else double(1),
		    if(mdata) jtmd else integer(1),
		    as.integer(ndyst),
		    as.integer(1),# jalg = 1 <==> AGNES
		    meth,# integer
		    integer(n),
		    ner = integer(n),
		    ban = double(n),
		    ac = as.double(0),
                    par.method,
		    merge = matrix(0:0, n - 1, 2), # integer
                    DUP = FALSE,
		    PACKAGE = "cluster")
    if(!diss) {
	##give warning if some dissimilarities are missing.
	if(res$ok == -1)
	    stop("No clustering performed, NA-values in the dissimilarity matrix.\n" )
        if(keep.diss) {
            ## adapt Fortran output to S:
            ## convert lower matrix,read by rows, to upper matrix, read by rows.
            disv <- res$dis[-1]
            disv[disv == -1] <- NA
            disv <- disv[upper.to.lower.tri.inds(n)]
            class(disv) <- dissiCl
            attr(disv, "Size") <- nrow(x)
            attr(disv, "Metric") <- metric
            attr(disv, "Labels") <- dimnames(x)[[1]]
        }
	##add labels to Fortran output
	if(length(dimnames(x)[[1]]) != 0)
	    order.lab <- dimnames(x)[[1]][res$ner]
    }
    else {
        if(keep.diss) disv <- x
	##add labels to Fortran output
	if(length(attr(x, "Labels")) != 0)
	    order.lab <- attr(x, "Labels")[res$ner]
    }
    clustering <- list(order = res$ner, height = res$ban[-1], ac = res$ac,
		       merge = res$merge, diss = if(keep.diss)disv,
                       call = match.call(), method = METHODS[meth])
    if(exists("order.lab"))
	clustering$order.lab <- order.lab
    if(keep.data && !diss) {
	if(mdata) x2[x2 == valmisdat] <- NA
	clustering$data <- x2
    }
    class(clustering) <- c("agnes", "twins")
    clustering
}

summary.agnes <- function(object, ...)
{
    class(object) <- "summary.agnes"
    object
}

print.agnes <- function(x, ...)
{
    cat("Call:	", deparse(x$call),
	"\nAgglomerative coefficient: ", format(x$ac, ...),
	"\nOrder of objects:\n")
    print(if(length(x$order.lab) != 0) x$order.lab else x$order,
	  quote = FALSE, ...)
    cat("Height (summary):\n");		print(summary(x$height), ...)
    cat("\nAvailable components:\n");	print(names(x), ...)
    invisible(x)
}

print.summary.agnes <- function(x, ...)
{
    ## a bit more than print.agnes() ..
    cat("Object of class `agnes' from call:\n", deparse(x$call),
	"\nAgglomerative coefficient: ", format(x$ac, ...),
	"\nOrder of objects:\n")
    print(if(length(x$order.lab) != 0) x$order.lab else x$order,
	  quote = FALSE, ...)
    cat("Merge:\n");			print(x$merge, ...)
    cat("Height:\n");			print(x$height, ...)
    if(!is.null(x$diss)) { ## Dissimilarities:
	cat("\n");			print(summary(x$diss, ...))
    }
    cat("\nAvailable components:\n");	print(names(x), ...)
    invisible(x)
}

### $Id: clara.q,v 1.23 2004/03/11 16:29:19 maechler Exp $

#### CLARA := Clustering LARge Applications
####
#### Note that the algorithm is O(n), but O(ns^2) where ns == sampsize

clara <- function(x, k, metric = "euclidean", stand = FALSE,
		  samples = 5, sampsize = min(n, 40 + 2 * k),
		  trace = 0, keep.data = TRUE, keepdata,
                  rngR = FALSE)
{
    ## check type of input matrix and values of input numbers
    x <- data.matrix(x)
    if(!is.numeric(x)) stop("x is not a numeric dataframe or matrix.")
    n <- nrow(x)
    if((k <- as.integer(k)) < 1 || k > n - 1)
	stop("The number of cluster should be at least 1 and at most n-1." )
    if((sampsize <- as.integer(sampsize)) < max(2,k))
	stop(paste(c("'sampsize' should be at least", max(2,k),
		     " = max(2, number of clusters)"), collapse = " "))
    if(n < sampsize)
	stop(paste(c("`sampsize' =", sampsize,
                     "should not be larger than the number of objects,", n),
		   collapse = " "))
    if((samples <- as.integer(samples)) < 1)
	stop("'samples' should be at least 1")

    ## `keepdata' is just for back-compatibility; "keep.*" is R-like
    if(!missing(keepdata) && missing(keep.data)) {
        warning("argument name `keepdata' is deprecated;",
                " use `keep.data' instead")
        keep.data <- keepdata
    }
    namx <- dimnames(x)[[1]]
    ## standardize, if necessary
    data <- x2 <- if(stand) scale(x, scale = apply(x, 2, meanabsdev)) else x
    ## put info about metric, size and NAs in arguments for the .C call
    jp <- ncol(x2)

    if((mdata <- any(inax <- is.na(x2)))) { # TRUE if x[] has any NAs
	jtmd <- as.integer(ifelse(apply(inax, 2, any), -1, 1))
	## VALue for MISsing DATa
	valmisdat <- 1.1* max(abs(range(x2, na.rm=TRUE)))
	x2[inax] <- valmisdat
    }

    if((trace <- as.integer(trace)))
	cat("calling .C(\"clara\", *):\n")
    res <- .C("clara",
	      n,
	      jp,
	      k,
	      clu = as.double(x2),
	      nran  = samples,
	      nsam  = sampsize,
	      dis   = double(1 + (sampsize * (sampsize - 1))/2),
	      mdata = as.integer(mdata),
	      valmd = if(mdata) rep(valmisdat, jp) else -1.,
	      jtmd  = if(mdata) jtmd else integer(1),
	      ndyst = as.integer(if(metric == "manhattan") 2 else 1),
              as.logical(rngR),
	      integer(sampsize),# = nrepr
	      integer(sampsize),# = nsel
	      sample= integer(sampsize),# = nbest
	      integer(k),		# = nr
	      med = integer(k),		# = nrx
	      double(k),		# = radus
	      double(k),		# = ttd
	      double(k),		# = ratt
	      avdis  = double(k),	# = ttbes
	      maxdis = double(k),	# = rdbes
	      ratdis = double(k),	# = rabes
	      size  = integer(k),	# = mtt
	      obj   = double(1),
	      avsil = double(k),
	      ttsil = double(1),
	      silinf = matrix(0, sampsize, 4),
	      jstop = integer(1),
	      trace = trace,
	      tmp  = double (3 * sampsize),
	      itmp = integer(6 * sampsize),
	      DUP = FALSE,
	      PACKAGE = "cluster")
    ## give a warning when errors occured
    if(res$jstop) {
	if(mdata && any(aNA <- apply(inax,1, all))) {
	    i <- which(aNA)
	    stop("Observation",
		 if(length(i) > 1)
		 paste("s  c(",paste(i, collapse=","),")  have", sep='')
		 else paste("", i, "has"),
		 " *only* NAs --> omit for clustering")
	} ## else
	if(res$jstop == 1)
	    stop("Each of the random samples contains objects between which\n",
		 " no distance can be computed.")
	if(res$jstop == 2)
	    stop("For each of the ", samples,
		 " samples, at least one object was found which\n could not",
		 " be assigned to a cluster (because of missing values).")
	## else {cannot happen}
	stop("invalid 'jstop' from .C(\"clara\",.): ", res$jstop)
    }
    sildim <- res$silinf[, 4]
    ## adapt C output to S:
    ## convert lower matrix, read by rows, to upper matrix, read by rows.
    disv <- res$dis[-1]
    disv[disv == -1] <- NA
    disv <- disv[upper.to.lower.tri.inds(sampsize)]
    class(disv) <- dissiCl
    attr(disv, "Size") <- sampsize
    attr(disv, "Metric") <- metric
    attr(disv, "Labels") <- namx[res$sample]
    ## add labels to C output
    res$med <- x[res$med, , drop = FALSE]
    res$clu <- as.integer(res$clu[1:n])
    if(!is.null(namx)) {
	sildim <- namx[sildim]
	res$sample <- namx[res$sample]
	names(res$clu) <- namx
    }
    ## add dimnames to C output
    r <- list(sample = res$sample, medoids = res$med,
	      clustering = res$clu, objective = res$obj,
	      clusinfo = cbind(size = res$size, "max_diss" = res$maxdis,
	      "av_diss" = res$avdis, isolation = res$ratdis),
	      diss = disv, call = match.call())
    if(k > 1) {
	dimnames(res$silinf) <- list(sildim,
				     c("cluster", "neighbor", "sil_width", ""))
	r$silinfo <- list(widths = res$silinf[, -4],
			  clus.avg.widths = res$avsil,
			  avg.width = res$ttsil)
    }
    if(keep.data) r$data <- data
    class(r) <- c("clara", "partition")
    r
}

print.clara <- function(x, ...)
{
    cat("Call:	", deparse(x$call),
	"\nMedoids:\n");		print(x$medoids, ...)
    cat("Objective function:\t ", format(x$objective, ...),"\n",
	"Clustering vector: \t", sep=""); str(x$clustering, vec.len = 7)
    cat("Cluster sizes:	    \t", x$clusinfo[,1],
	"\nBest sample:\n");		print(x$sample, quote = FALSE, ...)
    cat("\nAvailable components:\n");	print(names(x), ...)
    invisible(x)
}

summary.clara <- function(object, ...)
{
    class(object) <- "summary.clara"
    object
}

print.summary.clara <- function(x, ...)
{
    cat("Object of class `clara' from call:\n", deparse(x$call),
	"\nMedoids:\n");		print(x$medoids, ...)
    cat("Objective function:\t ", format(x$objective, ...),
	"\nNumerical information per cluster:\n")
    print(x$clusinfo, ...)
    if(has.sil <- !is.null(x$silinfo)) {
	cat("Average silhouette width per cluster:\n")
	print(x$silinfo[[2]], ...)
	cat("Average silhouette width of best sample:",
	    format(x$silinfo[[3]], ...), "\n")
    }
    cat("\nBest sample:\n");		print(x$sample, quote = FALSE, ...)
    cat("Clustering vector:\n");	print(x$clustering, ...)
    if(has.sil) {
	cat("\nSilhouette plot information for best sample:\n")
	print(x$silinfo[[1]], ...)
    }
    if(!is.null(x$diss)) { ## Dissimilarities:
	cat("\n");			print(summary(x$diss, ...))
    }
    cat("\nAvailable components:\n");	print(names(x), ...)
    invisible(x)
}

#### $Id: daisy.q,v 1.16 2004/03/11 16:26:40 maechler Exp maechler $
daisy <-
function(x, metric = c("euclidean","manhattan"), stand = FALSE, type = list())
{
    ## check type of input matrix
    if(length(dx <- dim(x)) != 2 || !(is.data.frame(x) || is.numeric(x)))
	stop("x is not a dataframe or a numeric matrix.")
    n <- dx[1]# nrow
    p <- dx[2]# ncol
    varnms <- dimnames(x)[[2]]
    pColl <- function(n) paste(n, collapse = ", ")
    if(length(type)) {
	if(!is.list(type) || is.null(ntyp <- names(type)) || any(ntyp == ""))
            stop("invalid ", sQuote("type"),"; must be named list")
        ## check each component to be valid column names or numbers:
        for(nt in ntyp) {
            cvec <- type[[nt]]
            if(is.character(cvec)) {
                if(!is.null(varnms) && !all(cvec %in% varnms))
                    stop("type$", nt, " has invalid column names")
            }
            else if(is.numeric(cvec)) {
                if(!all(1 <= cvec & cvec <= p))
                    stop("type$", nt, " must be in 1:ncol(x)")
            }
            else stop("type$", nt, " must contain column names or numbers")
        }
	tA <- type$asymm
	tS <- type$symm
	if((!is.null(tA) || !is.null(tS))) {
            ## tA and tS might be character and integer!
	    d.bin <- cbind(as.data.frame(x[, tA, drop= FALSE]),
                                         x[, tS, drop= FALSE])
            lenB <- sapply(lapply(d.bin, function(y)
                                 levels(as.factor(y))), length)
	    if(any(lenB > 2))
		stop("at least one binary variable has more than 2 levels.")
	    if(any(lenB < 2))
		warning("at least one binary variable has not 2 different levels.")
	    ## Convert factors to integer, such that ("0","1") --> (0,1):
	    if(any(is.f <- sapply(d.bin, is.factor)))
		d.bin[is.f] <- lapply(d.bin[is.f],
				      function(f) as.integer(as.character(f)))
	    if(!all(sapply(d.bin, function(y)
			   is.logical(y) ||
			   all(sort(unique(as.numeric(y[!is.na(y)])))%in% 0:1))))
		stop("at least one binary variable has values other than 0,1, and NA")
	}
    }
    ## transform variables and construct `type' vector
    if(is.data.frame(x)) {
	type2 <- sapply(x, data.class)
	x <- data.matrix(x)
    } else type2 <- rep("numeric", p)
    if(length(type)) {
	tT <- type$ ordratio
	tL <- type$ logratio
	x[, names(type2[tT])] <- unclass(as.ordered(x[, names(type2[tT])]))
	x[, names(type2[tL])] <- log10(		    x[, names(type2[tL])])
	type2[tA] <- "A"
	type2[tS] <- "S"
	type2[tT] <- "T" # was "O" (till 2000-12-14) accidentally !
    }
    type2[tI <- type2 %in% c("numeric", "integer") ] <- "I"
    if(any(tI) && any(iBin <- apply(x[,tI, drop = FALSE],2,
				    function(v) length(table(v)) == 2)))
	warning("binary variable(s) ", pColl(which(tI)[iBin]),
		" treated as interval scaled")

    type2[type2 == "ordered"] <- "O"
    type2[type2 == "factor"] <- "N"
    if(any(ilog <- type2 == "logical")) {
	warning("setting `logical' variable",if(sum(ilog)>1)"s " else " ",
		pColl(which(ilog)), " to type `asymm'")
	type2[ilog] <- "A"
    }
    ## standardize, if necessary
    if(all(type2 == "I")) {
	if(stand) {
            x <- scale(x, center = TRUE, scale = FALSE) #-> 0-means
            sx <- colMeans(abs(x))
            if(any(sx == 0)) {
                warning(sQuote("x"), " has constant columns ",
                        pColl(which(sx == 0)), "; these are standardized to 0")
                sx[sx == 0] <- 1
            }
	    x <- scale(x, center = FALSE, scale = sx)
        }
	jdat <- 2
	metric <- match.arg(metric)
	ndyst <- if(metric == "manhattan") 2 else 1
    }
    else { ## mixed case
	if(!missing(metric))
	    warning("`metric' is not used with mixed variables")
        colR <- apply(x, 2, range, na.rm = TRUE)
	colmin <- colR[1,]
	sx <- colR[2,] - colmin
        if(any(sx == 0))
            sx[sx == 0] <- 1
	x <- scale(x, center = colmin, scale = sx)
	jdat <- 1
	ndyst <- 0
    }
    ##	type2 <- paste(type2, collapse = "")
    typeCodes <- c('A','S','N','O','I','T')
    type3 <- match(type2, typeCodes)# integer
    if(any(ina <- is.na(type3)))
	stop("invalid type ", type2[ina],
             " for column numbers ", pColl(which(is.na)))
    if((mdata <- any(inax <- is.na(x)))) { # TRUE if x[] has any NAs
	jtmd <- as.integer(ifelse(apply(inax, 2, any), -1, 1))
	## VALue for MISsing DATa
	valmisdat <- 1.1* max(abs(range(x, na.rm=TRUE)))
	x[inax] <- valmisdat
	valmd <- rep(valmisdat, p)
    }
    ## call Fortran routine
    storage.mode(x) <- "double"
    disv <- .Fortran("daisy",
		     n,
		     p,
		     x,
		     if(mdata)valmd else double(1),
		     if(mdata) jtmd else integer(1),
		     as.integer(jdat),
		     type3,		# vtype
		     as.integer(ndyst),
		     dis = double((n * (n - 1))/2),
                     NAOK = TRUE,# only to allow "+- Inf"
		     DUP = FALSE,
		     PACKAGE = "cluster")$dis
    ## adapt Fortran output to S:
    ## convert lower matrix, read by rows, to upper matrix, read by rows.
    disv[disv == -1] <- NA
    full <- matrix(0, n, n)
    full[!lower.tri(full, diag = TRUE)] <- disv
    disv <- t(full)[lower.tri(full)]
    ## give warning if some dissimilarities are missimg
    if(any(is.na(disv))) attr(disv, "NA.message") <-
	"NA-values in the dissimilarity matrix !"
    ## construct S object -- "dist" methods are *there* !
    class(disv) <- dissiCl
    attr(disv, "Labels") <- dimnames(x)[[1]]
    attr(disv, "Size") <- n
    attr(disv, "Metric") <- ifelse(!ndyst, "mixed", metric)
    if(!ndyst) attr(disv, "Types") <- typeCodes[type3]
    disv
}

print.dissimilarity <- function(x, ...)
{
    cat("Dissimilarities :\n")
    print(as.vector(x), ...)
    cat("\n")
    if(!is.null(attr(x, "na.message")))
	cat("Warning : ", attr(x, "NA.message"), "\n")
    cat("Metric : ", attr(x, "Metric"),
	if(!is.null(aT <- attr(x,"Types")))
	paste(";  Types =", paste(aT, collapse=", ")), "\n")
    cat("Number of objects : ", attr(x, "Size"), "\n", sep="")
    invisible(x)
}

summary.dissimilarity <- function(object, ...)
{
    sx <- summary(as.vector(object), ...)
    at <- attributes(object)
    r <- c(list(summ = sx, n = length(object)), at[names(at) != "class"])
    class(r) <- "summary.dissimilarity"
    r
}

print.summary.dissimilarity <- function(x, ...)
{
    cat(x$n, "dissimilarities, summarized :\n")
    print(x$summ, ...)
    cat("Metric : ", x $ Metric,
	if(!is.null(aT <- x $ Types))
	paste(";  Types =", paste(aT, collapse=", ")), "\n")
    cat("Number of objects : ", x $ Size, "\n", sep="")
    if(!is.null(x $ na.message))
	cat("Warning : ", x $ NA.message, "\n")
    invisible(x)
}
### $Id: diana.q,v 1.16 2004/03/11 16:26:40 maechler Exp maechler $

diana <- function(x, diss = inherits(x, "dist"),
		  metric = "euclidean", stand = FALSE,
                  keep.diss = n < 100, keep.data = !diss)
{
   if((diss <- as.logical(diss))) {
	## check type of input vector
	if(any(is.na(x))) stop(..msg$error["NAdiss"])
	if(data.class(x) != "dissimilarity") { # try to convert to
	    if(!is.null(dim(x))) {
		x <- as.dist(x) # or give an error
	    } else {
		## possibly convert input *vector*
		if(!is.numeric(x) || is.na(n <- sizeDiss(x)))
		    stop(..msg$error["non.diss"])
		attr(x, "Size") <- n
	    }
	    class(x) <- dissiCl
	    if(is.null(attr(x,"Metric"))) attr(x, "Metric") <- "unspecified"
	}
	n <- as.integer(attr(x, "Size"))
	dv <- x[lower.to.upper.tri.inds(n)]
	## prepare arguments for the Fortran call
	dv <- c(0., dv)# double
	jp <- as.integer(1)
	mdata <- FALSE
	ndyst <- 0
	x2 <- double(1)
    }
    else {
	## check input matrix and standardize, if necessary
	x <- data.matrix(x)
	if(!is.numeric(x)) stop("x is not a numeric dataframe or matrix.")
	x2 <- if(stand) scale(x, scale = apply(x, 2, meanabsdev)) else x
	ndyst <- if(metric == "manhattan") 2 else 1
	n <- nrow(x2)
	jp <- ncol(x2)
	if((mdata <- any(inax <- is.na(x2)))) { # TRUE if x[] has any NAs
	    jtmd <- as.integer(ifelse(apply(inax, 2, any), -1, 1))
	    ## VALue for MISsing DATa
	    valmisdat <- 1.1* max(abs(range(x2, na.rm=TRUE)))
	    x2[inax] <- valmisdat
	    valmd <- rep(valmisdat, jp)
	}
	dv <- double(1 + (n * (n - 1))/2)
    }
    res <- .Fortran("twins",
		    n,
		    jp,
		    as.double(x2),
		    dv,
		    dis = double(if(keep.diss) length(dv) else 1),
		    ok = as.integer(if(keep.diss) diss + 10 else diss),
		    if(mdata)valmd else double(1),
		    if(mdata) jtmd else integer(1),
		    as.integer(ndyst),
		    as.integer(2),# jalg = 2 <==> DIANA
		    as.integer(0),# ~ method
		    integer(n),
		    ner = integer(n),
		    ban = double(n),
		    dc = as.double(0),
                    double(1),
		    merge = matrix(0:0, n - 1, 2), # integer
                    DUP = FALSE,
		    PACKAGE = "cluster")
    if(!diss) {
	## give warning if some dissimilarities are missing.
	if(res$ok == -1)
	    stop("No clustering performed, NA's in dissimilarity matrix.\n")
        if(keep.diss) {
            ## adapt Fortran output to S:
            ## convert lower matrix, read by rows, to upper matrix, read by rows.
            disv <- res$dis[-1]
            disv[disv == -1] <- NA
            disv <- disv[upper.to.lower.tri.inds(n)]
            class(disv) <- dissiCl
            attr(disv, "Size") <- nrow(x)
            attr(disv, "Metric") <- metric
            attr(disv, "Labels") <- dimnames(x)[[1]]
        }
	## add labels to Fortran output
	if(length(dimnames(x)[[1]]) != 0)
	    order.lab <- dimnames(x)[[1]][res$ner]
    }
    else {
        if(keep.diss) disv <- x
	## add labels to Fortran output
	if(length(attr(x, "Labels")) != 0)
	    order.lab <- attr(x, "Labels")[res$ner]
    }
    clustering <- list(order = res$ner, height = res$ban[-1], dc = res$dc,
		       merge = res$merge, diss = if(keep.diss)disv,
                       call = match.call())
    if(exists("order.lab"))
	clustering$order.lab <- order.lab
    if(keep.data && !diss) {
	if(mdata) x2[x2 == valmisdat] <- NA
	clustering$data <- x2
    }
    class(clustering) <- c("diana", "twins")
    clustering
}

print.diana <- function(x, ...)
{
    cat("Merge:\n")
    print(x$merge, ...)
    cat("Order of objects:\n")
    print(if (length(x$order.lab) != 0) x$order.lab else x$order,
	  quote = FALSE, ...)
    cat("Height:\n")
    print(x$height, ...)
    cat("Divisive coefficient:\n")
    print(x$dc, ...)
    cat("\nAvailable components:\n")
    print(names(x), ...)
    invisible(x)
}

summary.diana <- function(object, ...)
{
    class(object) <- "summary.diana"
    object
}

print.summary.diana <- function(x, ...)
{
    cat("Merge:\n");			print(x$merge, ...)
    cat("Order of objects:\n")
    print(if(length(x$order.lab)) x$order.lab else x$order, quote = FALSE, ...)
    cat("Height:\n");			print(x$height, ...)
    cat("Divisive coefficient:\n");	print(x$dc, ...)
    if(!is.null(x$diss)) { ## Dissimilarities:
	cat("\n");			print(summary(x$diss, ...))
    }
    cat("\nAvailable components:\n");	print(names(x), ...)
    invisible(x)
}
####-*- mode: R; kept-old-versions: 12;  kept-new-versions: 20; -*-

#### ellipsoidhull : Find (and optionally draw)
#### -----------   the smallest ellipsoid containining a set of points
####
#### Just making the algorithms in clusplot() available more generally
#### ( --> ./plotpart.q )

### Author: Martin Maechler, Date: 21 Jan 2002, 15:41

ellipsoidhull <-
    function(x, tol = 0.01, maxit = 5000,
             ret.wt = FALSE, ret.sqdist = FALSE, ret.pr = FALSE)
{
    if(!is.matrix(x) || !is.numeric(x))
        stop("`x' must be numeric  n x p matrix")
    if(any(is.na(x))) {
        warning("omitting NAs")
        x <- na.omit(x)
    }
    n <- nrow(x)
    if(n == 0) stop("no points without missing values")
    p <- ncol(x)

    res <- .C("spannel",
              n,
              ndep= p,
              dat = cbind(1., x),
              sqdist = double(n),
              l1 = double((p+1) ^ 2),
              double(p),
              double(p),
              prob = double(n),
              double(p+1),
              eps = as.double(tol),
              maxit = as.integer(maxit),
              ierr = integer(1),# 0 or non-zero
              PACKAGE = "cluster")
    if(res$ierr != 0)
        cat("Error in Fortran routine computing the spanning ellipsoid,",
            "\n probably collinear data\n", sep="")
    if(any(res$prob < 0) || all(res$prob == 0))
        stop("computed some negative or all 0 `prob'abilities")
    conv <- res$maxit  < maxit
    if(!conv)
        warning("possibly not converged in ", maxit, " iterations")
    conv <- conv && res$ierr == 0

    cov <- cov.wt(x, res$prob)
    ## cov.wt() in R has extra wt[] scaling; revert here
    res <- list(loc = cov$center,
                cov = cov$cov * (1 - sum(cov$wt^2)),
                d2  = weighted.mean(res$sqdist, res$prob),
                wt  = if(ret.wt) cov$wt,
                sqdist = if(ret.sqdist) res$sqdist,
                prob= if(ret.pr) res$prob,
                tol = tol,
                eps = max(res$sqdist) - p,
                it  = res$maxit,
                maxit= maxit,
                ierr = res$ierr,
                conv = conv)

    class(res) <- "ellipsoid"
    res
}

print.ellipsoid <- function(x, digits = max(1, getOption("digits") - 2), ...)
{
    d <- length(x$loc)
    cat("`ellipsoid' in", d, "dimensions:\n center = (",
        format(x$loc, digits=digits),
        "); squared ave.radius d^2 = ", format(x$d2, digits=digits),
        "\n and shape matrix =\n")
    print(x$cov, digits = digits, ...)
    cat("  hence,",if(d==2)"area" else "volume"," = ",
        format(volume(x), digits=digits),"\n")
    if(!x$conv) {
        cat("\n** Warning: ** the algorithm did not terminate reliably!\n  ",
            if(x$ierr) "most probably because of collinear data"
            else "(in the available number of iterations)", "\n")
    }
    invisible(x)
}

volume <- function(object) UseMethod("volume")
volume.ellipsoid <- function(object) {
    A <- object$cov
    pi * object$d2 * sqrt(det(A))
}

## For p = 2 :
##   Return (x[i],y[i]) points, i = 1:n, on boundary of ellipse, given
##   by 2 x 2 matrix A[], origin `loc' and d(xy, loc) ^2 = `d2'
ellipsoidPoints <- function(A, d2, loc, n = 201)
{
    if(length(d <- dim(A)) != 2 || (p <- d[1]) != d[2])
        stop("`A' must be p x p  cov-matrix defining an ellipsoid")
    if(p == 2) {
        detA <- A[1, 1] * A[2, 2] - A[1, 2]^2
        yl2 <- A[2, 2] * d2
        y <- seq( - sqrt(yl2), sqrt(yl2), leng = n)
        sqrt.discr <- sqrt(detA * pmax(0, yl2 - y^2))/A[2, 2]
        sqrt.discr[c(1, n)] <- 0
        b <- loc[1] + A[1, 2]/A[2, 2] * y
        y <- loc[2] + y
        return(rbind(cbind(    b - sqrt.discr,      y),
                     cbind(rev(b + sqrt.discr), rev(y))))
    } else { ## p >= 3
        detA <- det(A)
        ##-- need something like polar coordinates
        stop("ellipsoidPoints() not yet implemented for p >= 3 dim.")
    }
}

predict.ellipsoid <- function(object, n.out = 201, ...)
    ellipsoidPoints(object$cov, d2 = object$d2, loc= object$loc, n = n.out)
#### $Id: fanny.q,v 1.15 2004/03/11 16:26:40 maechler Exp maechler $
fanny <- function(x, k, diss = inherits(x, "dist"),
		  metric = "euclidean", stand = FALSE)
{
   if((diss <- as.logical(diss))) {
	## check type of input vector
	if(any(is.na(x))) stop(..msg$error["NAdiss"])
	if(data.class(x) != "dissimilarity") { # try to convert to
	    if(!is.null(dim(x))) {
		x <- as.dist(x) # or give an error
	    } else {
		## possibly convert input *vector*
		if(!is.numeric(x) || is.na(n <- sizeDiss(x)))
		    stop(..msg$error["non.diss"])
		attr(x, "Size") <- n
	    }
	    class(x) <- dissiCl
	    if(is.null(attr(x,"Metric"))) attr(x, "Metric") <- "unspecified"
	}
	## prepare arguments for the Fortran call
	n <- attr(x, "Size")
	dv <- as.double(c(x, 0))
	jp <- 1
	mdata <- FALSE
	ndyst <- 0
	x2 <- double(n)
	jdyss <- 1
    }
    else {
	## check input matrix and standardize, if necessary
	x <- data.matrix(x)
	if(!is.numeric(x)) stop("x is not a numeric dataframe or matrix.")
	x2 <- if(stand) scale(x, scale = apply(x, 2, meanabsdev)) else x
	## put info about metric, size and NAs in arguments for the Fortran call
	ndyst <- if(metric == "manhattan") 2 else 1
	n <- nrow(x2)
	jp <- ncol(x2)
	if((mdata <- any(inax <- is.na(x2)))) { # TRUE if x[] has any NAs
	    jtmd <- as.integer(ifelse(apply(inax, 2, any), -1, 1))
	    ## VALue for MISsing DATa
	    valmisdat <- 1.1* max(abs(range(x2, na.rm=TRUE)))
	    x2[inax] <- valmisdat
	    valmd <- rep(valmisdat, jp)
	}
	dv <- double(1 + (n * (n - 1))/2)
	jdyss <- 0
    }
    if((k <- as.integer(k)) < 1 || k > n%/%2 - 1)
	stop("`k' (number of clusters) must be in {1,2, .., n/2 -1}")
    ## call Fortran routine
    storage.mode(x2) <- "double"
    res <- .Fortran("fanny",
		    as.integer(n),
		    as.integer(jp),
		    k,
		    x2,
		    dis = dv,
		    ok = as.integer(jdyss),
		    if(mdata)valmd else double(1),
		    if(mdata) jtmd else integer(1),
		    as.integer(ndyst),
		    integer(n),
		    integer(n),
		    integer(n),
		    double(n),
		    p = matrix(0., n, k),
		    matrix(0., n, k),
		    avsil = double(k),
		    integer(k),
		    double(k),
		    double(k),
		    double(n),
		    ttsil = as.double(0),
		    eda = as.double(0),
		    edb = as.double(0),
		    obj = double(2),
		    clu = integer(n),
		    silinf = matrix(0., n, 4),
		    as.double(1e-15),
		    PACKAGE = "cluster")
    sildim <- res$silinf[, 4]
    if(diss) {
	disv <- x
	## add labels to Fortran output
	if(length(attr(x, "Labels")) != 0) {
	    sildim <- attr(x, "Labels")[sildim]
	    dimnames(res$p) <- list(attr(x, "Labels"), NULL)
	    names(res$clu) <- attr(x, "Labels")
	}
    }
    else {
	## give warning if some dissimilarities are missing.
	if(res$ok == -1)
	    stop("No clustering performed, NA-values in the dissimilarity matrix.")
	disv <- res$dis[ - (1 + (n * (n - 1))/2)]
	disv[disv == -1] <- NA
	class(disv) <- dissiCl
	attr(disv, "Size") <- nrow(x)
	attr(disv, "Metric") <- metric
	attr(disv, "Labels") <- dimnames(x)[[1]]
	## add labels to Fortran output
	if(length(dimnames(x)[[1]]) != 0) {
	    sildim <- dimnames(x)[[1]][sildim]
	    dimnames(res$p) <- list(dimnames(x)[[1]], NULL)
	    names(res$clu) <- dimnames(x)[[1]]
	}
    }
    ## add dimnames to Fortran output
    names(res$obj) <- c("iterations", "objective")
    res$coeff <- c(res$eda, res$edb)
    names(res$coeff) <- c("dunn_coeff", "normalized")

    r <- list(membership = res$p, coeff = res$coeff,
		       clustering = res$clu, objective = res$obj,
		       diss = disv, call = match.call())
    if(k != 1) {
	dimnames(res$silinf) <- list(sildim,
				     c("cluster", "neighbor", "sil_width", ""))
	r$silinfo <- list(widths = res$silinf[, -4],
                          clus.avg.widths = res$avsil[1:k],
                          avg.width = res$ttsil)
    }
    if(!diss) {
	if(mdata) x2[x2 == valmisdat] <- NA
	r$data <- x2
    }
    class(r) <- c("fanny", "partition")
    r
}

print.fanny <- function(x, ...)
{
    print(x$objective, ...)
    cat("Membership coefficients:\n")
    print(x$membership, ...)
    cat("Coefficients:\n")
    print(x$coeff, ...)
    cat("Closest hard clustering:\n")
    print(x$clustering, ...)
    cat("\nAvailable components:\n")
    print(names(x), ...)
    invisible(x)
}

summary.fanny <- function(object, ...)
{
    class(object) <- "summary.fanny"
    object
}

print.summary.fanny <- function(x, ...)
{
    print(x$objective, ...)
    cat("Membership coefficients:\n");	print(x$membership, ...)
    cat("Coefficients:\n");		print(x$coeff, ...)
    cat("Closest hard clustering:\n");	print(x$clustering, ...)
    if(length(x$silinfo) != 0) {
	cat("\nSilhouette plot information:\n")
	print(x$silinfo[[1]], ...)
	cat("Average silhouette width per cluster:\n")
	print(x$silinfo[[2]], ...)
	cat("Average silhouette width of total data set:\n")
	print(x$silinfo[[3]], ...)
    }
    if(!is.null(x$diss)) { ## Dissimilarities:
	cat("\n");			print(summary(x$diss, ...))
    }
    cat("\nAvailable components:\n");	print(names(x), ...)
    invisible(x)
}
#### Cluster - Internal Utilities
#### ============================ (new by Martin Mchler)

## This was size(); seems slightly useful in general
sizeDiss <- function(d)
{
    ## find `n' for d == dissimilarity-like(<n obs.>), i.e. length(d)= n(n-1)/2
    discr <- 1 + 8 * length(d)
    sqrtdiscr <- round(sqrt(discr))
    if(sqrtdiscr^2 == discr) (1 + sqrtdiscr)/2 else NA
}

## used in  ./agnes.q, ./clara.q,  ./diana.q  und ./pam.q :

lower.to.upper.tri.inds <- function(n)
{
    n1 <- as.integer(n - 1)
    if(n1 < 1) stop("`n' must be >= 2")
    else if(n1 == 1) 1:1
    else rep(1:n1, 1:n1) +
        c(0, unlist(lapply(2:n1, function(k) cumsum(c(0, (n - 2):(n - k))))))
}

upper.to.lower.tri.inds <- function(n)
{
    if((n2 <- as.integer(n - 2)) < 0) stop("`n' must be >= 2")
    rep(1 + cumsum(0:n2), (n - 1):1) +
	unlist(lapply(0:n2, function(k) cumsum(k:n2)))
}

#### consider to *not* export these when I will use a name space :

meanabsdev <- function(y) mean(abs(y - mean(y, na.rm = TRUE)), na.rm = TRUE)
####-*- mode: R; kept-old-versions: 12;  kept-new-versions: 20; -*-

mona <- function(x)
{
    levs <- function(y) levels(as.factor(y))

    ## check type of input matrix
    if(!is.matrix(x) && !is.data.frame(x))
        stop("x must be a matrix or data frame.")
    if(!all(sapply(lapply(as.data.frame(x), levs), length) == 2))
        stop("All variables must be binary (factor with 2 levels).")
    n <- nrow(x)
    jp <- ncol(x)
    ## change levels of input matrix
    x2 <- apply(as.matrix(x), 2, factor)
    x2[x2 == "1"] <- "0"
    x2[x2 == "2"] <- "1"
    x2[is.na(x2)] <- "2"
    ##	x2 <- paste(x2, collapse = "")
    ##	storage.mode(x2) <- "character"
    ## call Fortran routine
    storage.mode(x2) <- "integer"
    res <- .Fortran("mona",
                    as.integer(n),
                    as.integer(jp),
                    x2 = x2,# x[,]
                    error = as.integer(0),
                    nban = integer(n),
                    ner = integer(n),
                    integer(n),
                    lava = integer(n),
                    integer(jp),
                    PACKAGE = "cluster")
    ## stop with a message when two many missing values:
    if(res$error != 0) {
        ch <- "No clustering performed, "
        switch(res$error,
               ## 1 :
               stop(ch,"an object was found with all values missing."),
               ## 2 :
               stop(ch,"a variable was found with at least 50% missing values."),
               ## 3 :
               stop(ch,"a variable was found with all non missing values identical."),
               ## 4 :
               stop(ch,"all variables have at least one missing value.")
               )
    }
    ##O res$x2 <- matrix(as.numeric(substring(res$x2,
    ##O                                      1:nchar(res$x2), 1:nchar(res$x2))),
    ##O                      n, jp)
    storage.mode(res$x2) <- "integer" # keeping dim()
    dimnames(res$x2) <- dnx <- dimnames(x)
    ## add labels to Fortran output
    if(length(dnx[[2]]) != 0) {
        lava <- as.character(res$lava)
        lava[lava != "0"] <- dnx[[2]][res$lava]
        lava[lava == "0"] <- "NULL"
        res$lava <- lava
    }
    ## construct "mona" object
    clustering <- list(data = res$x2, order = res$ner,
                       variable = res$lava[ -1 ], step = res$nban[-1],
                       call = match.call())
    if(length(dnx[[1]]) != 0)
        clustering$order.lab <- dnx[[1]][res$ner]
    class(clustering) <- "mona"
    clustering
}

print.mona <- function(x, ...)
{
    cat("Revised data:\n")
    print(x$data, quote = FALSE, ...)
    cat("Order of objects:\n")
    print(if (length(x$order.lab) != 0) x$order.lab else x$order,
          quote = FALSE, ...)
    cat("Variable used:\n")
    print(x$variable, quote = FALSE, ...)
    cat("Separation step:\n")
    print(x$step, ...)
    cat("\nAvailable components:\n")
    print(names(x), ...)
    invisible(x)
}

## FIXME: These should differ from print()

summary.mona <- function(object, ...)
{
    class(object) <- "summary.mona"
    object
}

print.summary.mona <- function(x, ...)
{
    print.mona(x, ...)
    invisible(x)
}

#### PAM : Partitioning Around Medoids
#### --- $Id: pam.q,v 1.19 2004/03/11 16:26:40 maechler Exp maechler $
pam <- function(x, k, diss = inherits(x, "dist"),
		metric = "euclidean", stand = FALSE, cluster.only = FALSE,
                keep.diss = !diss && !cluster.only && n < 100,
                keep.data = !diss && !cluster.only)
{
   if((diss <- as.logical(diss))) {
	## check type of input vector
	if(any(is.na(x))) stop(..msg$error["NAdiss"])
	if(data.class(x) != "dissimilarity") { # try to convert to
	    if(!is.null(dim(x))) {
		x <- as.dist(x) # or give an error
	    } else {
		## possibly convert input *vector*
		if(!is.numeric(x) || is.na(n <- sizeDiss(x)))
		    stop(..msg$error["non.diss"])
		attr(x, "Size") <- n
	    }
	    class(x) <- dissiCl
	    if(is.null(attr(x,"Metric"))) attr(x, "Metric") <- "unspecified"
	}
	## adapt S dissimilarities to Fortran:
	## convert upper matrix, read by rows, to lower matrix, read by rows.
	n <- attr(x, "Size")
	dv <- x[lower.to.upper.tri.inds(n)]
	## prepare arguments for the Fortran call
	dv <- c(0, dv)
	jp <- 1
	mdata <- FALSE
	ndyst <- 0
	x2 <- double(n)
    }
    else {
	## check input matrix and standardize, if necessary
	x <- data.matrix(x)
	if(!is.numeric(x)) stop("x is not a numeric dataframe or matrix.")
	x2 <- if(stand) scale(x, scale = apply(x, 2, meanabsdev)) else x
	## put info about metric, size and NAs in arguments for the Fortran call
	ndyst <- if(metric == "manhattan") 2 else 1
	n <- nrow(x2)
	jp <- ncol(x2)
	if((mdata <- any(inax <- is.na(x2)))) { # TRUE if x[] has any NAs
	    jtmd <- as.integer(ifelse(apply(inax, 2, any), -1, 1))
	    ## VALue for MISsing DATa
	    valmisdat <- 1.1* max(abs(range(x2, na.rm=TRUE)))
	    x2[inax] <- valmisdat
	    valmd <- rep(valmisdat, jp)
	}
	dv <- double(1 + (n * (n - 1))/2)
    }
    if((k <- as.integer(k)) < 1 || k >= n)
	stop("Number of clusters `k' must be in {1,2, .., n-1}; hence n >= 2")
    ## call Fortran routine
    storage.mode(dv) <- "double"
    storage.mode(x2) <- "double"
    res <- .C("pam",
	      as.integer(n),
	      as.integer(jp),
	      k,
	      x = x2,
	      dys = dv,
	      jdyss = as.integer(diss),
	      if(mdata)valmd else double(1),
	      if(mdata) jtmd else integer(1),
	      as.integer(ndyst),
	      integer(n),		# nsend[]
	      logical(n),		# nrepr[]
	      integer(if(cluster.only) 1 else n), # nelem[]
	      double(n),		# radus[]
	      double(n),		# damer[]
	      avsil = double(n),	# `ttd'
	      double(n),		# separ[]
	      ttsil = as.double(0),
	      obj = c(as.double(cluster.only), 0.),# in & out!
	      med = integer(if(cluster.only) 1 else k),
	      clu = integer(n),
	      clusinf = if(cluster.only) 0. else matrix(0., k, 5),
	      silinf  = if(cluster.only) 0. else matrix(0., n, 4),
	      isol = integer(if(cluster.only) 1 else k),
	      DUP = FALSE, # care!!
	      PACKAGE = "cluster")

    xLab <- if(diss) attr(x, "Labels") else dimnames(x)[[1]]
    if(length(xLab) > 0)
        names(res$clu) <- xLab
    if(cluster.only)
        return(res$clu)

    ## Else, usually
    sildim <- res$silinf[, 4]
    if(diss) {
	disv <- x
	## add labels to Fortran output
	if(length(xLab) > 0) {
	    sildim <- xLab[sildim]
	    res$med <- xLab[res$med]
	}
    }
    else {
	## give warning if some dissimilarities are missing.
	if(res$jdyss == -1)
	    stop("No clustering performed, NAs in the computed dissimilarity matrix.")
        if(keep.diss) {
            ## adapt Fortran output to S:
            ## convert lower matrix, read by rows, to upper matrix, read by rows.
            disv <- res$dys[-1]
            disv[disv == -1] <- NA
            disv <- disv[upper.to.lower.tri.inds(n)]
            class(disv) <- dissiCl
            attr(disv, "Size") <- nrow(x)
            attr(disv, "Metric") <- metric
            attr(disv, "Labels") <- dimnames(x)[[1]]
        }
	## add labels to Fortran output
	res$med <- x[res$med,  , drop =FALSE]
	if(length(xLab) > 0)
	    sildim <- xLab[sildim]
    }
    ## add dimnames to Fortran output
    names(res$obj) <- c("build", "swap")
    res$isol <- factor(res$isol, levels = 0:2, labels = c("no", "L", "L*"))
    names(res$isol) <- 1:k
    dimnames(res$clusinf) <- list(NULL, c("size", "max_diss", "av_diss",
					  "diameter", "separation"))
    ## construct S object
    r <-
	list(medoids = res$med, clustering = res$clu,
	     objective = res$obj, isolation = res$isol,
	     clusinfo = res$clusinf,
	     silinfo = if(k != 1) {
		 dimnames(res$silinf) <-
		     list(sildim, c("cluster", "neighbor", "sil_width", ""))
		 list(widths = res$silinf[, -4],
		      clus.avg.widths = res$avsil[1:k],
		      avg.width = res$ttsil)
	     },
	     diss = if(keep.diss)disv,
	     call = match.call())
    if(keep.data && !diss) {
	if(mdata) x2[x2 == valmisdat] <- NA
	r$data <- x2
    }
    class(r) <- c("pam", "partition")
    r
}

print.pam <- function(x, ...)
{
    cat("Medoids:\n")
    print(x$medoids, ...)
    cat("Clustering vector:\n")
    print(x$clustering, ...)
    cat("Objective function:\n")
    print(x$objective, ...)
    cat("\nAvailable components:\n")
    print(names(x), ...)
    invisible(x)
}

summary.pam <- function(object, ...)
{
    class(object) <- "summary.pam"
    object
}

print.summary.pam <- function(x, ...)
{
    cat("Medoids:\n");			print(x$medoids, ...)
    cat("Clustering vector:\n");	print(x$clustering, ...)
    cat("Objective function:\n");	print(x$objective, ...)
    cat("\nNumerical information per cluster:\n"); print(x$clusinfo, ...)
    cat("\nIsolated clusters:\n L-clusters: ")
    print(names(x$isolation[x$isolation == "L"]), quote = FALSE, ...)
    cat(" L*-clusters: ")
    print(names(x$isolation[x$isolation == "L*"]), quote = FALSE, ...)
    if(length(x$silinfo) != 0) {
	cat("\nSilhouette plot information:\n")
	print(x$silinfo[[1]], ...)
	cat("Average silhouette width per cluster:\n")
	print(x$silinfo[[2]], ...)
	cat("Average silhouette width of total data set:\n")
	print(x$silinfo[[3]], ...)
    }
    if(!is.null(x$diss)) { ## Dissimilarities:
	cat("\n");			print(summary(x$diss, ...))
    }
    cat("\nAvailable components:\n");	print(names(x), ...)
    invisible(x)
}

### $Id: plothier.q,v 1.14 2003/04/30 14:12:11 maechler Exp $

pltree <- function(x, ...) UseMethod("pltree")

## note: pltree() can have an `xlab' in "..." (plot.hclust has an explicit one)
pltree.twins <- function(x, main = paste("Dendrogram of ", deparse(x$call)),
			 labels = NULL, ylab = "Height", ...)
{
    if(is.null(labels) && length(x$order.lab) != 0)
	labels <- x$order.lab[sort.list(x$order)]

    ## calling plot.hclust() via generic :
    plot(structure(list(merge = x$merge, order = x$order,
                        height = sort(x$height), labels = labels,
                        call = x$call, method = x$method),
                   class = "hclust"),
         main = main, ylab = ylab, ...)
    invisible()
}

bannerplot <-
function(x, w = rev(x$height), fromLeft = TRUE,
	 main, sub, xlab = "Height", adj = 0, col = c(2, 0), border = 0,
	 axes = TRUE, frame.plot = axes, rev.xax = !fromLeft, xax.pretty = TRUE,
	 labels = NULL, nmax.lab = 35, max.strlen = 5,
	 yax.do = axes && length(x$order) <= nmax.lab,
	 yaxRight = fromLeft, y.mar = 2.4 + max.strlen / 2.5, ...)
{
    m <- max(w)
    if(axes) {
	if(xax.pretty) {
	    n <- if(!is.logical(xax.pretty)) xax.pretty else formals(pretty)$n
	    at.vals <- pretty(c(0,w), n = n)
	    n <- length(at.vals <- at.vals[at.vals <= m])
	    if(at.vals[n] * 1.01 < m) {
		lab.vals <- c(at.vals, signif(m, 3))
		at.vals <-  c(at.vals, m)
	    } else lab.vals <- at.vals
	} else { # old default for plot.agnes() and plot.diana()
	    ss <- seq(0, floor(m), length = 11)# => intervals = 1/10 {total}
	    at.vals  <- c(ss, m)
	    lab.vals <- round(at.vals, 2)
	}
    }
    if(fromLeft) {
	w <- rbind(w, m - w)
	if(missing(col)) col <- rev(col)
    } else { ## from Right
	w <- rbind(m - w, w)
	if(axes && rev.xax) {
	    at.vals <- m - rev(at.vals)## == c(0, ss + m - floor(m))
	    lab.vals <- rev(lab.vals)
	}
    }
    if(yax.do) {
	ax <- if(yaxRight)
	    list(side = 4, pos = m)
	else list(side = 2, pos = 0)
	if((pm <- par("mar"))[ax$side] < y.mar) {
	    ## need space besides y axis for labeling
	    pm[ax$side] <- y.mar
	    op <- par(mar = pm)
	    on.exit(par(op))
	}
    }
    barplot(w, xlab = xlab, horiz = TRUE, space = 0, axes = FALSE,
	    col = col, border = border, mgp = c(2.5, 1, 0), ...)
    if(frame.plot && (border == 0 || border == par("bg")))
	rect(0, 0, m, ncol(w))

    title(main = main, sub = sub, adj = adj)
    if(axes) {
	axis(1, at = at.vals, labels = lab.vals, ...)
	if(yax.do) {
	    if(is.null(labels))
		labels <- rev(if (length(x$order.lab) != 0)
			      substring(x$order.lab, 1,max.strlen) else x$order)
	    axis(ax$side, at = 0:(length(x$order) - 1), las = 1,
		 labels = labels, pos = ax$pos, mgp = c(3, 1.25, 0), ...)
	}
    }
}

## plot.diana() [further down] & plot.agnes() are  almost identical;
## --  made  bannerplot() a stand-alone function
## --> maybe *merge* these two into one plot.twins()
plot.agnes <-
function(x, ask = FALSE, which.plots = NULL, main = NULL,
	 sub = paste("Agglomerative Coefficient = ", round(x$ac, digits = 2)),
         adj = 0, nmax.lab = 35, max.strlen = 5, xax.pretty = TRUE, ...)
{
    if(is.null(main)) {
	## Different default for banner & pltree:
	cl <- deparse(x$call)
	main1 <- paste("Banner of ", cl)
	main2 <- paste("Dendrogram of ", cl)
    }
    else { # same title for both
	main1 <- main2 <- main
    }

    if(is.null(which.plots) && !ask)
	which.plots <- 1:2
    if(ask && is.null(which.plots)) { ## Use `menu' ..
	tmenu <- paste("plot ", ## choices :
		       c("All", "Banner", "Clustering Tree"))
	do.all <- FALSE
	repeat {
	    if(!do.all)
		pick <- menu(tmenu, title =
			     "\nMake a plot selection (or 0 to exit):\n") + 1
	    switch(pick,
		   return(invisible()), # 0 -> exit loop
		   do.all <- TRUE,# 1 : All
		   bannerplot(x, fromLeft = TRUE,
			      main = main1, sub = sub, adj = adj,
			      xax.pretty = 10,
			      nmax.lab= nmax.lab, max.strlen= max.strlen, ...),
		   pltree (x, main = main2, sub = sub, ...) # 3
		   )
	    if(do.all) { pick <- pick + 1; do.all <- pick <= length(tmenu) + 1}
	}
    }
    else {
	ask <- prod(par("mfcol")) < length(which.plots) && dev.interactive()
	if(ask) {
	    op <- par(ask = TRUE)
	    on.exit(par(op))
	}
	for(i in which.plots)
	switch(i,
	       bannerplot(x, fromLeft = TRUE,
			  main = main1, sub = sub, adj = adj,
			  xax.pretty = 10,
			  nmax.lab = nmax.lab, max.strlen = max.strlen, ...),
	       pltree (x, main = main2, sub = sub, ...)
	       )
    }
    invisible()
}

plot.diana <-
function(x, ask = FALSE, which.plots = NULL, main = NULL,
	 sub  = paste("Divisive Coefficient = ", round(x$dc, digits = 2)),
	 adj = 0, nmax.lab = 35, max.strlen = 5, xax.pretty = TRUE, ...)
{
    if(is.null(main)) {
	## Different default for banner & pltree:
	cl <- deparse(x$call)
	main1 <- paste("Banner of ", cl)
	main2 <- paste("Dendrogram of ", cl)
    }
    else { # same title for both
	main1 <- main2 <- main
    }

    if(is.null(which.plots) && !ask)
	which.plots <- 1:2
    if(ask && is.null(which.plots)) { ## Use `menu' ..
	tmenu <- paste("plot ", ## choices :
		       c("All", "Banner", "Clustering Tree"))
	do.all <- FALSE
	repeat {
	    if(!do.all)
		pick <- menu(tmenu, title =
			     "\nMake a plot selection (or 0 to exit):\n") + 1
	    switch(pick,
		   return(invisible()), # 0 -> exit loop
		   do.all <- TRUE,# 1 : All
		   bannerplot(x, fromLeft = FALSE,
			      main = main1, sub = sub, adj = adj,
			      xax.pretty = 10,
			      nmax.lab= nmax.lab, max.strlen= max.strlen, ...),
		   pltree (x, main = main2, sub = sub, ...)
		   )
	    if(do.all) { pick <- pick + 1; do.all <- pick <= length(tmenu) + 1}
	}
    }
    else {
	ask <- prod(par("mfcol")) < length(which.plots) && dev.interactive()
	if(ask) {
	    op <- par(ask = TRUE)
	    on.exit(par(op))
	}
	for(i in which.plots)
	switch(i,
	       bannerplot(x, fromLeft = FALSE, main = main1, sub = sub,
			  adj = adj, xax.pretty = 10,
			  nmax.lab = nmax.lab, max.strlen = max.strlen, ...),# 1
	       pltree (x, main = main2, sub = sub, ...) # i = 2
	       )
    }
    invisible()
}

plot.mona <- function(x, main = paste("Banner of ", deparse(x$call)),
		      sub = NULL, xlab = "Separation step",
		      col = c(2,0), axes = TRUE, adj = 0,
		      nmax.lab = 35, max.strlen = 5, ...)
{
    w <- rev(x$step)
    m <- max(w)
    if(any(i0 <- w == 0))
	w[i0] <- m <- m+1
    bannerplot(x[c("order","order.lab")], w = w, fromLeft = TRUE,
	       yaxRight = FALSE, col = col, main = main, sub = sub, xlab = xlab,
	       adj= adj, axes= axes, nmax.lab= nmax.lab, max.strlen= max.strlen,
	       xax.pretty = m+1, ...)
    names <- paste(" ", rev(x$variable))
    is.na(names) <- i0
    text(w, 1:length(names) - 0.5, names, adj = 0, col = col[1], ...)
}
### $Id: plotpart.q,v 1.22 2004/03/08 10:56:34 maechler Exp maechler $
plot.partition <-
function(x, ask = FALSE, which.plots = NULL,
         nmax.lab = 40, max.strlen = 5, data = x$data, dist = NULL,
         cor = TRUE, stand = FALSE, lines = 2,
         shade = FALSE, color = FALSE, labels = 0, plotchar = TRUE,
         span = TRUE, xlim = NULL, ylim = NULL, main = NULL, ...)
{
    if(is.null(x$data))# data not kept
        x$data <- data
    if(is.null(x$data) && !is.null(dist))
        x$diss <- dist
    if(is.null(which.plots) && !ask)
        which.plots <- {
            if(is.null(x$data) && (is.null(x$diss) || inherits(x, "clara")))
                2 ## no clusplot
            else 1:2
        }
    if(ask && is.null(which.plots)) { ## Use `menu' ..
        tmenu <- paste("plot ", ## choices :
                       c("All", "Clusplot", "Silhouette Plot"))
        do.all <- FALSE
        repeat {
            if(!do.all)
                pick <- menu(tmenu, title =
                             "\nMake a plot selection (or 0 to exit):\n") + 1
            switch(pick,
                   return(invisible())# 0 -> exit loop
                   ,
                   do.all <- TRUE# 1 : All
                   ,
                   clusplot(x, cor = cor, stand = stand, lines = lines,
                            shade = shade, color = color, labels = labels,
                            plotchar = plotchar, span = span,
                            xlim = xlim, ylim = ylim, main = main, ...)
                   ,
                   plot(silhouette(x), nmax.lab, max.strlen, main = main)
                   )
            if(do.all) { pick <- pick + 1; do.all <- pick <= length(tmenu) + 1}
        }
    }
    else {
        ask <- prod(par("mfcol")) < length(which.plots) && dev.interactive()
        if(ask) { op <- par(ask = TRUE); on.exit(par(op)) }
        for(i in which.plots)
        switch(i,
               clusplot(x, cor = cor, stand = stand, lines = lines,
                        shade = shade, color = color, labels = labels,
                        plotchar = plotchar, span = span,
                        xlim = xlim, ylim = ylim, main = main, ...)
               ,
               plot(silhouette(x), nmax.lab, max.strlen, main = main)
               )
    }
    invisible()
}

clusplot <- function(x, ...) UseMethod("clusplot")

clusplot.default <-
function(x, clus, diss = FALSE, cor = TRUE, stand = FALSE, lines = 2,
         shade = FALSE, color = FALSE, labels = 0, plotchar = TRUE,
         col.p = "dark green", # was 5 (= shaded col)
         col.txt = col.p, col.clus = if(color) c(2, 4, 6, 3) else 5,
         span = TRUE, xlim = NULL, ylim = NULL,
         main = paste("CLUSPLOT(", deparse(substitute(x)),")"),
         verbose = getOption("verbose"),
         ...)
{
    force(main)
    if(is.data.frame(x))
        x <- data.matrix(x)
    if(!is.numeric(x))
        stop("x is not numeric")

    if(diss) {
        if(any(is.na(x)))
            stop("NA-values are not allowed in dist-like `x'.")
        if(inherits(x, "dist")) {
            n <- attr(x, "Size")
            labels1 <- attr(x, "Labels")
        }
        else { # x (num.vector or square matrix) must be transformed into diss.
            siz <- sizeDiss(x)
            if(is.na(siz)) {
                if((n <- nrow(x)) != ncol(x))
                    stop("Distances must be result of dist or a square matrix.")
                if(all.equal(x, t(x)) != TRUE)
                    stop("the square matrix is not symmetric.")
                labels1 <- dimnames(x)[[1]]
            }
            else {
                if(!is.vector(x)) {
                    labels1 <- attr(x, "Labels") # possibly NULL
                    x <- as.matrix(x)
                    if((n <- nrow(x)) == ncol(x) && all.equal(x, t(x)) == TRUE) {
                        labels1 <- dimnames(x)[[1]]
                    }
                    else {
                        ## Hmm, when does this ever happen :
                        ## numeric, not-dist, non-vector, not symmetric matrix ?
                        warning(">>>>> funny case in clusplot.default() -- please report!\n")
                        ## if(n != sizeDiss(x)) ...
                        if(is.null(labels1))
                            labels1 <- 1:sizeDiss(x)
                        attr(x, "Size") <- sizeDiss(x)
                    }
                }
                else {
                    attr(x, "Size") <- n <- siz
                    labels1 <- 1:n
                }
            }
        }
        if(is.null(labels1)) labels1 <- 1:n

        if(paste(R.version$major, R.version$minor, sep=".") < 1.5) {
            ## a simplified (add = T) version of R 1.5's cmdscale()
            cmdscale <- function (d, k = 2, add = TRUE, ...) {
                if (any(is.na(d)))
                    stop("NA values not allowed in d")
                if (is.null(n <- attr(d, "Size"))) {
                    d <- as.matrix(d)
                    x <- d^2
                    if ((n <- nrow(x)) != ncol(x))
                        stop("Distances must be result of dist or a square matrix")
                }
                else {
                    x <- matrix(0, n, n)
                    if(add) d0 <- x
                    x[row(x) > col(x)] <- d^2
                    x <- x + t(x)
                    if(add) {
                        d0[row(x) > col(x)] <- d
                        d <- d0 + t(d0)
                    }
                }
                storage.mode(x) <- "double"
                x <- .C("dblcen", x=x, as.integer(n), PACKAGE="mva")$x
                if(add) { ## solve the additive constant problem
                    i2 <- n + (i <- 1:n)
                    Z <- matrix(0, 2*n, 2*n)
                    Z[cbind(i2,i)] <- -1
                    Z[ i, i2] <- -x
                    Z[i2, i2] <- .C("dblcen", x= 2*d, as.integer(n),PACKAGE="mva")$x
                    e <- eigen(Z,symmetric = FALSE, only.val = TRUE)$values
                    add.c <- max(Re(e))
                    x <- matrix(double(n*n), n, n)
                    non.diag <- row(d) != col(d)
                    x[non.diag] <- (d[non.diag] + add.c)^2
                }
                e <- eigen(-x/2, symmetric = TRUE)
                ev <- e$values[1:k]
                points <- e$vectors[, 1:k] %*% diag(sqrt(ev), k)
                rn <- if(is.matrix(d)) rownames(d) else names(d)
                dimnames(points) <- list(rn, NULL)
                evalus <- e$values[-n]
                list(points = points, eig = ev, ac = if(add) add.c else 0,
                     GOF = sum(ev)/c(sum(abs(evalus)),
                     sum(evalus[evalus > 0])))
            }
        } ## cmdscale() -- if R version < 1.5

        x1 <- cmdscale(x, k = 2, eig = TRUE, add = TRUE)
        if(x1$ac < 0)
            x1 <- cmdscale(x, k = 2, eig = TRUE)
        var.dec <- x1$GOF[2] # always in [0,1]
        x1 <- x1$points
    }
    else { ## Not (diss)
        if(!is.matrix(x)) stop("x is not a data matrix")
        if(any(is.na(x))) {
            y <- is.na(x)
            if(any(apply(y, 1, all)))
                stop("one or more objects contain only missing values")
            if(any(apply(y, 2, all)))
                stop("one or more variables contain only missing values")
            x <- apply(x, 2, function(x)
                   { x[is.na(x)] <- median(x, na.rm = TRUE); x } )
            cat("Missing values were displaced by the median of the corresponding variable(s)\n")

        }

        n <- nrow(x)
        labels1 <- dimnames(x)[[1]]
        if(is.null(labels1)) labels1 <- 1:n

        if(ncol(x) == 1) {
            hulp <- rep(0, length(x))
            x1 <- matrix(c(t(x), hulp), ncol = 2)
            var.dec <- 1
        }
        else {
            prim.pr <- princomp(x, scores = TRUE, cor = ncol(x) != 2)
            var.dec <- cumsum(prim.pr$sdev^2/sum(prim.pr$ sdev^2))[2]
            x1 <- prim.pr$scores
            x1 <- cbind(x1[, 1], x1[, 2])
        }
    }

    ## --- The 2D space is setup and points are in x1[,]  (aantal x 2) ---

    clus <- as.vector(clus)
    if(length(clus) != n)
        stop("The clustering vector has not the good length")
    clus <- as.factor(clus)
    if(any(is.na(clus)))
        stop("NA-values are not allowed in clustering vector")
    if(stand)
        x1 <- scale(x1)

    rangx <- range(x1[, 1])
    rangy <- range(x1[, 2])
    minx <- rangx[1]
    maxx <- rangx[2]
    miny <- rangy[1]
    maxy <- rangy[2]
    levclus <- levels(clus)
    nC <- length(levclus) # the number of clusters
    z <- A <- vector("list", nC)
    maxima <- loc <- matrix(0, nrow = nC, ncol = 2)
    d2 <- verhoud <- numeric(nC)
    verhouding <- 0
    ## num1 .. num6 : all used only once -- there are more constants anyway
    num3 <- 90
    num6 <- 70

    for(i in 1:nC) { ##-------------  i-th cluster  --------------
	x <- x1[clus == levclus[i],, drop = FALSE ]
        aantal <- nrow(x) # number of observations in cluster [i]
        cov <- var(if(aantal == 1) {
                     if(verbose)
                         cat("cluster",i," has only one observation ..\n")
                     rbind(x, c(0, 0))
                   } else x)
        x.1 <- range(x[, 1])
        y.1 <- range(x[, 2])
        notrank2 <- qr(cov, tol = 0.001)$rank != 2
        if(!span && notrank2) {
            d2[i] <- 1
            if((abs(diff(x.1)) > (diff(rangx)/70)) ||
               (abs(diff(y.1)) > (diff(rangy)/50))) {
                loc[i, ] <- c(x.1[1] + diff(x.1)/2, y.1[1] + diff(y.1)/2)
                a <- sqrt((loc[i, 1] - x.1[1])^2 +
                          (loc[i, 2] - y.1[1])^2)
                a <- a + 0.05 * a
                num2 <- 40
                if(abs(diff(x.1)) > diff(rangx)/70 ) {
                    ind1 <- which.max(x[,1])
                    ind2 <- which.min(x[,1])
                    q <- atan((x[ind1, 2] - x[ind2, 2])/
                              (x[ind1, 1] - x[ind2, 1]))
                    b <-
                        if(diff(rangy) == 0)
                            1
                        else if(abs(diff(y.1)) > diff(rangy)/50)
                            diff(y.1)/10 ## num1 <- 10
                        else diff(rangy)/num2
                }
                else {
                    b <- if(diff(rangx) == 0) 1 else diff(rangx)/num2
                    q <- pi/2
                }
                D <- diag(c(a^2, b^2))
                R <- rbind(c(cos(q), -sin(q)),
                           c(sin(q),  cos(q)))
                A[[i]] <- (R %*% D) %*% t(R)
            }
            else {
                a <- diff(rangx)/num3
                b <- diff(rangy)/num6
                if(a == 0) a <- 1
                if(b == 0) b <- 1
                A[[i]] <- diag(c(a^2, b^2))
                loc[i, ] <- x[1, ]
            }
            oppervlak <- pi * a * b
        }
        else if(span && notrank2) {
            d2[i] <- 1
            if(sum(x[, 1] != x[1, 1]) != 0 ||
               sum(x[, 2] != x[1, 2]) != 0) {
                loc[i, ] <- c(x.1[1] + diff(x.1)/2,
                              y.1[1] + diff(y.1)/2)
                a <- sqrt((loc[i, 1] - x.1[1])^2 +
                          (loc[i, 2] - y.1[1])^2)
                if(any(x[, 1] != x[1, 1])) {
                    ind1 <- which.max(x[,1])
                    ind2 <- which.min(x[,1])
                    q <- atan((x[ind1, 2] - x[ind2, 2])/
                              (x[ind1, 1] - x[ind2, 1]))
                }
                else {
                    q <- pi/2
                }
                b <- 1e-7
                D <- diag(c(a^2, b^2))
                R <- rbind(c(cos(q), -sin(q)),
                           c(sin(q),  cos(q)))
                A[[i]] <- (R %*% D) %*% t(R)
            }
            else {
                a <- diff(rangx)/num3
                b <- diff(rangy)/num6
                if(a == 0) a <- 1
                if(b == 0) b <- 1
                A[[i]] <- diag(c(a^2, b^2))
                loc[i, ] <- x[1, ]
            }
            oppervlak <- pi * a * b

        }
        else { ## rank2
            if(!span) {
                loc[i, ] <- apply(x, 2, mean)
                d2[i] <- max(mahalanobis(x, loc[i, ], cov))
                ## * (1+ 0.01)^2  --- dropped factor for back-compatibility
            }
            else { ## span and rank2
                if(verbose)
                    cat("span & rank2 : calling \"spannel\" ..\n")
                k <- as.integer(2)
		res <- .C("spannel",
			  aantal,
			  ndep= k,
			  dat = cbind(1., x),
			  sqdist = double(aantal),
			  l1 = double((k+1) ^ 2),
			  double(k),
			  double(k),
			  prob = double(aantal),
			  double(k+1),
			  eps = (0.01),## convergence tol.
			  maxit = as.integer(5000),
			  ierr = integer(1),
			  PACKAGE = "cluster")
                if(res$ierr != 0)
                    ## MM : exactmve not available here !
                    cat("Error in Fortran routine for the spanning ellipsoid,",
                        "\n rank problem??\n", sep="")

                cov <- cov.wt(x, res$prob)
                loc[i, ] <- cov$center
                ## NB: cov.wt() in R has extra wt[] scaling; revert here:
                cov <- cov$cov * (1 - sum(cov$wt^2))
                d2[i] <- weighted.mean(res$sqdist, res$prob)

                if(verbose)
                    cat("ellipse( A= (", format(cov[1,]),"*", format(cov[2,2]),
                        "),\n\td2=", format(d2[i]),
                        ", loc[]=", format(loc[i, ]), ")\n")
            }
            A[[i]] <- cov
            ## oppervlak (flam.)  =  area (Engl.)
            oppervlak <- pi * d2[i] * sqrt(cov[1, 1] * cov[2, 2] - cov[1, 2]^2)
        }

        z[[i]] <- ellipsoidPoints(A[[i]], d2[i], loc[i, ], n= 201)
        maxima[i, ] <- z[[i]][201, ]
        rx <- range(z[[i]][, 1])
        ry <- range(z[[i]][, 2])
        minx <- min(minx, rx[1])
        maxx <- max(maxx, rx[2])
        miny <- min(miny, ry[1])
        maxy <- max(maxy, ry[2])
        verhoud[i] <- aantal/oppervlak
        if(verhoud[i] < 1e7)
            verhouding <- verhouding + verhoud[i]
    } ## end for( i-th cluster )

    if(verhouding == 0)
        verhouding <- 1
    ## num4 <- 37 ; num5 <- 3 --- but `41' is another constant
    density <- 3 + (verhoud * 37)/verhouding
    density[density > 41] <- 41
    if (span) {
        if (rangx[1] == rangx[2]) { ## diff(rangx)==0 : x-coords all the same
            minx <- x1[1, 1] - 1
            maxx <- x1[1, 1] + 1
        }
        if (rangy[1] == rangy[2]) { ## diff(rangy)==0 : y-coords all the same
            miny <- x1[1, 2] - 1
            maxy <- x1[1, 2] + 1
        }
    }
    if(is.null(xlim)) xlim <- c(minx,maxx)
    if(is.null(ylim)) ylim <- c(miny,maxy)
    if(length(col.p) < n) col.p <- rep(col.p, length= n)

    ## --- Now plotting starts ---

    plot(x1[, 1], x1[, 2], xlim = xlim, ylim = ylim,
         xlab = "Component 1", ylab = "Component 2",
         main = main,
         type = if(plotchar) "n" else "p", # if(plotchar) add points later
         col = col.p, ...)
    title(sub = paste("These two components explain",
          round(100 * var.dec, digits = 2), "% of the point variability."),
          adj = 0)

    if(color) {
        if(length(col.clus) < min(4,nC))
            stop("`col.clus' should have length 4 when color is TRUE")
        i.verh <- order(verhoud)
        jInd <- if(nC > 4) pam(verhoud[i.verh], 4)$clustering else 1:nC
        for(i in 1:nC) {
            k <- i.verh[i]
            polygon(z[[k]], density = if(shade) density[k] else 0,
                    col = col.clus[jInd[i]], ...)
        }
        col.clus <- col.clus[jInd][order(i.verh)]
    }
    else {
        for(i in 1:nC)
            polygon(z[[i]], density = if(shade) density[i] else 0,
                    col = col.clus, ...)
    }

    ## points after polygon in order to write ON TOP:
    if(plotchar) {
        karakter <- 1:19
        for(i in 1:nC) {
            iC <- clus == levclus[i]
            x <- x1[iC, , drop = FALSE]
            il <- 1+(i-1) %% 19
            points(x[, 1], x[, 2], pch = karakter[il], col = col.p[iC], ...)
        }
    }

    if((lines == 1 || lines == 2) && nC > 1) {
        ## Draw lines between all pairs of the  nC  cluster (centers)

        ## utilities for computing ellipse intersections:
        clas.snijpunt <- function(x, loc, m, n, p)
        {
            if(     loc[n, m] <= x[1, m] && x[1, m] <= loc[p, m]) x[1, ]
            else if(loc[n, m] <= x[2, m] && x[2, m] <= loc[p, m]) x[2, ]
            else NA
        }
        coord.snijp1 <- function(x, gemid)
            x[2, 2] - 2 * x[1, 2] * gemid + x[1, 1] * gemid^2
        coord.snijp2 <- function(x, d2, y)
            ((x[1, 1] * x[2, 2] - x[1, 2]^2) * d2)/y
        coord.snijp3 <- function(xx, y, gemid)
        {
            sy <- sqrt(y)
            sy <- c(sy, -sy)
            cbind(xx[1] + sy,
                  xx[2] + gemid*sy)
        }

        afstand <- matrix(0, ncol = nC, nrow = nC)
        for(i in 1:(nC - 1)) {
            for(j in (i + 1):nC) {
                gemid <- (loc[j, 2] - loc[i, 2])/(loc[j, 1] - loc[i, 1])
                s0 <- coord.snijp1(A[[i]], gemid)
                b0 <- coord.snijp2(A[[i]], d2[i], s0)
                snijp.1 <- coord.snijp3(loc[i,], y=b0, gemid)
                s1 <- coord.snijp1(A[[j]], gemid)
                b1 <- coord.snijp2(A[[j]], d2[j], s1)
                snijp.2 <- coord.snijp3(loc[j,], y=b1, gemid)
                if(loc[i, 1] != loc[j, 1]) {
                    if(loc[i, 1] < loc[j, 1]) {
                        punt.1 <- clas.snijpunt(snijp.1, loc, 1, i, j)
                        punt.2 <- clas.snijpunt(snijp.2, loc, 1, i, j)
                    }
                    else {
                        punt.1 <- clas.snijpunt(snijp.1, loc, 1, j, i)
                        punt.2 <- clas.snijpunt(snijp.2, loc, 1, j, i)
                    }
                }
                else {
                    if(loc[i, 2] < loc[j, 2]) {
                        punt.1 <- clas.snijpunt(snijp.1, loc, 2, i, j)
                        punt.2 <- clas.snijpunt(snijp.2, loc, 2, i, j)
                    }
                    else {
                        punt.1 <- clas.snijpunt(snijp.1, loc, 2, j, i)
                        punt.2 <- clas.snijpunt(snijp.2, loc, 2, j, i)
                    }
                }
                if(is.na(punt.1[1]) || is.na(punt.2[1]) ||
                   (sqrt((punt.1[1] - loc[i, 1])^2 +
                         (punt.1[2] - loc[i, 2])^2) +
                    sqrt((punt.2[1] - loc[j, 1])^2 +
                         (punt.2[2] - loc[j, 2])^2)) >
                   sqrt((loc[j, 1] - loc[i, 1])^2 +
                        (loc[j, 2] - loc[i, 2])^2))
                {
                    afstand[i, j] <- NA
                }
                else if(lines == 1) {
                    afstand[i, j] <- sqrt((loc[i, 1] - loc[j, 1])^2 +
                                          (loc[i, 2] - loc[j, 2])^2)
                    segments(loc[i, 1], loc[i, 2],
                             loc[j, 1], loc[j, 2], col = 6, ...)
                }
                else { ## lines == 2
                    afstand[i, j] <- sqrt((punt.1[1] - punt.2[1])^2 +
                                          (punt.1[2] - punt.2[2])^2)
                    segments(punt.1[1], punt.1[2],
                             punt.2[1], punt.2[2], col = 6, ...)
                }
            }
        }
        afstand <- t(afstand) + afstand
    }
    else afstand <- NULL

    if(labels) {
        if(labels == 1) {
            for(i in 1:nC) { ## add cluster border points
                m <- nrow(z[[i]])
                ni <- length(ii <- seq(1, m, by = max(1, m %/% 40)))
                x1 <- rbind(x1, z[[i]][ii, ])
                labels1 <- c(labels1, rep(levclus[i], ni))
                ## identify() only allows one color:
                ##col.txt <- c(col.txt, rep(col.clus[if(color) i else 1], ni))
            }
            identify(x1, labels = labels1, col = col.txt[1])
        }
        else {
            Stext <- function(xy, labs, ...) {
                ## FIXME: these displacements are not quite ok!
                xy[, 1] <- xy[, 1] + (maxx - minx)/130
                xy[, 2] <- xy[, 2] + (maxy - miny)/50
                text(xy, labels = labs, ...)
            }
            if(labels == 3 || labels == 2)
                Stext(x1, labels1, col = col.txt)
            if(labels %in% c(2,4,5))
                Stext(maxima, levclus, font = 4, col = col.clus)
            if(labels == 5)
                identify(x1, labels = labels1, col = col.txt[1])
        }
    }
    density[density == 41] <- NA
    invisible(list(Distances = afstand, Shading = density))
}

clusplot.partition <- function(x, main = NULL, dist = NULL, ...)
{
    if(is.null(main) && !is.null(x$call))
	main <- paste("clusplot(",format(x$call),")", sep="")
    if(length(x$data) != 0 &&
       (!any(is.na(x$data)) || data.class(x) == "clara"))
	clusplot.default(x$data, x$clustering, diss = FALSE, main = main, ...)
    else if(!is.null(dist))
        clusplot.default(dist, x$clustering, diss = TRUE, main = main, ...)
    else if(!is.null(x$diss))
        clusplot.default(x$diss, x$clustering, diss = TRUE, main = main, ...)
    else { ## try to find "x$diss" by looking at the pam() call:
        if(!is.null(x$call)) {
            xD <- try(eval(x$call[[2]], envir = parent.frame()))
            if(inherits(xD, "try-error") || !inherits(xD, "dist"))
                stop("no diss nor data found, nor the original argument of ",
                     deparse(x$call))
            ## else
            ## warning("both `x$diss' and `dist' are empty; ",
            ##         "trying to find the first argument of ", deparse(x$call))
            clusplot.default(xD, x$clustering, diss = TRUE, main = main, ...)
        }
        else stop("no diss nor data found for clusplot()'")
    }
}
silhouette <- function(x, ...) UseMethod("silhouette")

## Accessor and more:
silhouette.partition <- function(x, ...) {
    r <- x$silinfo$widths
    if(is.null(r))
        stop("invalid partition object")
    attr(r, "Ordered") <- TRUE # (cluster <increasing>, s.i <decreasing>)
    attr(r, "call") <- x$call
    class(r) <- "silhouette"
    r
}

silhouette.default <- function(x, dist, dmatrix, ...) {
    cll <- match.call()
    if(!is.null(cl <- x$clustering)) x <- cl
    n <- length(x)
    if(!all(x == round(x))) stop("`x' must only have integer codes")
    k <- length(clid <- sort(unique(x)))
    if(k <= 1 || k >= n)
        return(NA)
    ## check dist/dmatrix
    if(missing(dist)) {
        if(missing(dmatrix))
            stop("Need either a dissimilarity `dist' or diss.matrix `dmatrix'")
        if(is.null(dm <- dim(dmatrix)) || length(dm) != 2 || !all(n == dm))
            stop("`dmatrix' is not a dissimilarity matrix compatible to `x'")
    } else { # `dist'
        dist <- as.dist(dist) # hopefully
        if(n != attr(dist, "Size"))
            stop("clustering `x' and dissimilarity `dist' are incompatible")
        dmatrix <- as.matrix(dist)# so we can apply(.) below
    }
    wds <- matrix(NA, n,3, dimnames =
                  list(names(x), c("cluster","neighbor","sil_width")))
    for(j in 1:k) { # j-th cluster:
        Nj <- sum(iC <- x == clid[j])
        wds[iC, "cluster"] <- clid[j]
        a.i <- if(Nj > 1) colSums(dmatrix[iC, iC])/(Nj - 1) else 0
                                        # length(a.i)= Nj
        ## minimal distances to points in all other clusters:
        diC <- rbind(apply(dmatrix[!iC, iC, drop = FALSE], 2,
                           function(r) tapply(r, x[!iC], mean)))# (k-1) x Nj
        minC <- max.col(-t(diC))
        wds[iC,"neighbor"] <- clid[-j][minC]
        b.i <- diC[cbind(minC, seq(along = minC))]
        s.i <- ifelse(a.i != b.i, (b.i - a.i) / pmax(b.i, a.i), 0)
        wds[iC,"sil_width"] <- s.i
    }
    attr(wds, "Ordered") <- FALSE
    attr(wds, "call") <- cll
    class(wds) <- "silhouette"
    wds
}


sortSilhouette <- function(object, ...) {
    if(attr(object,"Ordered")) return(object)
    ## Else :
    if(is.null(n <- nrow(object)) || n < 1)
        stop("invalid silhouette structure")
    if(is.null(rownames(object)))
        rownames(object) <- as.character(1:n)
    ## k <- length(clid <- sort(unique(cl <- object[,"cluster"])))# cluster ID s
    cl <- object[,"cluster"]
    r <- object[iOrd <- order(cl, - object[,"sil_width"]) , , drop = FALSE]
    attributes(r) <- attributes(object) # but:
    attr(r,"iOrd") <- iOrd # the ordering
    attr(r,"Ordered") <- TRUE
    r
}

summary.silhouette <- function(object, FUN = mean, ...)
{
    if(ncol(object) != 3) stop("invalid `silhouette' object")
    cl <- object[, "cluster"]
    si <- object[, "sil_width"]
    r <- list(si.summary = summary(si, ...),
              clus.avg.widths = tapply(si, cl, FUN),
              clus.sizes = table(cl),
              avg.width = FUN(si),
              call = attr(object,"call"),
              Ordered = attr(object,"Ordered"))
    class(r) <- "summary.silhouette"
    r
}

print.summary.silhouette <- function(x, ...)
{
    k <- length(csiz <- x$clus.sizes)
    cat("Silhouette of", sum(csiz), "units in", k, "clusters",
        if(!is.null(x$call)) paste("from", deparse(x$call)),
        ":\nCluster sizes and average silhouette widths:\n")
    cwid <- x$clus.avg.widths
    names(cwid) <- csiz
    print(cwid, ...)
    cat("Individual silhouette widths:\n")
    print(x$si.summary, ...)
    invisible(x)
}


## This was the internal function silhouPlot() in plot.partition() :
plot.silhouette <-
    function(x, nmax.lab = 40, max.strlen = 5,
	     main = NULL, sub = NULL,
	     xlab = expression("Silhouette width " * s[i]),
	     col = "gray", do.col.sort = length(col) > 1,
	     border = 0, cex.names = par("cex.axis"),
	     do.n.k = TRUE, do.clus.stat = TRUE, ...)
{
    if(!is.matrix(x) || ncol(x) != 3)
	stop("No valid silhouette information (#{clusters} =? 1)")
    n <- nrow(x)
    x <- sortSilhouette(x)
    s <- rev(x[, "sil_width"])
    space <- c(0, rev(diff(cli <- x[, "cluster"])))
    space[space != 0] <- 0.5 # gap between clusters
    names <- if(n < nmax.lab)
	substring(rev(rownames(x)), 1, max.strlen)
    if(is.null(main)) {
	main <- "Silhouette plot"
	if(!is.null(cll <- attr(x,"call"))) { # drop initial "silhouette":
	    if(!is.na(charmatch("silhouette", deparse(cll[[1]]))))
		cll[[1]] <- as.name("FF")
	    main <-  paste(main, "of", sub("^FF","", deparse(cll)))
	}
    }
    smry <- summary(x)
    k <- length(nj <- smry$clus.sizes) # k clusters
    if(is.null(sub))
	sub <- paste("Average silhouette width : ",
		     round(smry$avg.width, digits = 2))
    if(do.col.sort && (lc <- length(col)) > 1) {
	if(lc == k)# cluster wise coloring
	    col <- col[cli]
	else if(lc != n)
	    col <- rep(col, length = n)
	col <- rev(col[attr(x, "iOrd")])
    }
    y <- barplot(s, space = space, names = names, xlab = xlab,
		 xlim = c(min(0, min(s)), 1),
		 horiz = TRUE, las = 1, mgp = c(2.5, 1, 0),
		 col = col, border = border, cex.names = cex.names, ...)
    title(main = main, sub = sub, adj = 0)
    if(do.n.k) {
	mtext(paste("n =", n),	adj = 0)
	mtext(substitute(k ~~ "clusters" ~~ C[j], list(k=k)), adj= 1)
    }
    if(do.clus.stat) {
	mtext(expression(paste(j," :  ", n[j]," | ", ave[i %in% Cj] ~~ s[i])),
	      adj = 1.04, line = -1.2)
	y <- rev(y)
	for(j in 1:k) {
	    yj <- mean(y[cli == j])
	    text(1, yj, paste(j,":  ", nj[j]," | ",
			      format(smry$clus.avg.widths[j], digits = 2)),
		 xpd = NA, adj = 0.8)
	}
    }
}
## no S4 methodology here; speedup :
.noGenerics <- TRUE
