/******************************************************************************
* Triv_aux.c - auxiliary routine to interface to tri-variate rep.	      *
*******************************************************************************
* (C) Gershon Elber, Technion, Israel Institute of Technology                 *
*******************************************************************************
* Written by Gershon Elber, Sep. 94.					      *
******************************************************************************/

#include "triv_loc.h"

/*****************************************************************************
* DESCRIPTION:                                                               M
* Given a tri-variate, returns its parametric domain.			     M
*                                                                            *
* PARAMETERS:                                                                M
*   TV:           Trivariate function to consider.                           M
*   UMin, UMax:   U Domain of TV will be placed herein.                      M
*   VMin, VMax:   V Domain of TV will be placed herein.                      M
*   WMin, WMax:   W Domain of TV will be placed herein.                      M
*                                                                            *
* RETURN VALUE:                                                              M
*   void                                                                     M
*                                                                            *
* KEYWORDS:                                                                  M
*   TrivTVDomain, trivariates                                                M
*****************************************************************************/
void TrivTVDomain(TrivTVStruct *TV,
		  CagdRType *UMin,
		  CagdRType *UMax,
		  CagdRType *VMin,
		  CagdRType *VMax,
		  CagdRType *WMin,
		  CagdRType *WMax)
{
    int UOrder,	VOrder, WOrder, ULen, VLen, WLen;

    switch (TV -> GType) {
	case TRIV_TVBEZIER_TYPE:
	    *UMin = 0.0;
	    *UMax = 1.0;
	    *VMin = 0.0;
	    *VMax = 1.0;
	    *WMin = 0.0;
	    *WMax = 1.0;
	    break;
	case TRIV_TVBSPLINE_TYPE:
	    UOrder = TV -> UOrder;
	    VOrder = TV -> VOrder;
	    WOrder = TV -> WOrder;
	    ULen = TV -> ULength;
	    VLen = TV -> VLength;
	    WLen = TV -> WLength;

	    *UMin = TV -> UKnotVector[UOrder - 1];
	    *UMax = TV -> UKnotVector[ULen];
	    *VMin = TV -> VKnotVector[VOrder - 1];
	    *VMax = TV -> VKnotVector[VLen];
	    *WMin = TV -> WKnotVector[WOrder - 1];
	    *WMax = TV -> WKnotVector[WLen];
	    break;
	default:
	    TRIV_FATAL_ERROR(TRIV_ERR_UNDEF_GEOM);
	    break;
    }
}

/*****************************************************************************
* DESCRIPTION:                                                               M
* Given a tri-variate and a domain - validate it.			     M
*                                                                            *
* PARAMETERS:                                                                M
*   TV:       To make sure t is in its Dir domain.                           M
*   t:        Parameter value to verify.                                     M
*   Dir:      Direction. Either U or V or W.                                 M
*                                                                            *
* RETURN VALUE:                                                              M
*   CagdBType:    TRUE if in domain, FALSE otherwise.                        M
*                                                                            *
* KEYWORDS:                                                                  M
*   TrivParamInDomain, trivariates                                           M
*****************************************************************************/
CagdBType TrivParamInDomain(TrivTVStruct *TV, CagdRType t, TrivTVDirType Dir)
{
    CagdRType UMin, UMax, VMin, VMax, WMin, WMax;

    TrivTVDomain(TV, &UMin, &UMax, &VMin, &VMax, &WMin, &WMax);

    switch (Dir) {
	case TRIV_CONST_U_DIR:
	    return t >= UMin && t <= UMax;
	case TRIV_CONST_V_DIR:
	    return t >= VMin && t <= VMax;
	case TRIV_CONST_W_DIR:
	    return t >= WMin && t <= WMax;
	default:
	    TRIV_FATAL_ERROR(TRIV_ERR_WRONG_DOMAIN);
	    return FALSE;
    }
}

/*****************************************************************************
* DESCRIPTION:                                                               M
* Given a tri-variate and a domain - validate it.			     M
*                                                                            *
* PARAMETERS:                                                                M
*   TV:       To make sure (u, v, w) is in its domain.                       M
*   u, v, w:  To verify if it is in TV's parametric domain.                  M
*                                                                            *
* RETURN VALUE:                                                              M
*   CagdBType:   TRUE if in domain, FALSE otherwise.                         M
*                                                                            *
* KEYWORDS:                                                                  M
*   TrivParamsInDomain, trivariates                                          M
*****************************************************************************/
CagdBType TrivParamsInDomain(TrivTVStruct *TV,
			     CagdRType u,
			     CagdRType v,
			     CagdRType w)
{
    CagdRType UMin, UMax, VMin, VMax, WMin, WMax;

    TrivTVDomain(TV, &UMin, &UMax, &VMin, &VMax, &WMin, &WMax);

    return u >= UMin && u <= UMax &&
           v >= VMin && v <= VMax &&
	   w >= WMin && w <= WMax;
}

/*****************************************************************************
* DESCRIPTION:                                                               M
* Given a tri-variate, returns a sub-region of it.                           M
*                                                                            *
* PARAMETERS:                                                                M
*   TV:        To extract asub-region from.                                  M
*   t1, t1:    Domain to extract from TV, in parametric direction Dir.       M
*   Dir:       Direction to extract the sub-region. Either U or V or W.      M
*                                                                            *
* RETURN VALUE:                                                              M
*   TrivTVStruct *:   A sub-region of TV from t1 to t2 in direction Dir.     M
*                                                                            *
* KEYWORDS:                                                                  M
*   TrivTVRegionFromTV, trivariates                                          M
*****************************************************************************/
TrivTVStruct *TrivTVRegionFromTV(TrivTVStruct *TV,
				 CagdRType t1,
				 CagdRType t2,
				 TrivTVDirType Dir)
{
    CagdRType RMin, RMax, SMin, SMax, TMin, TMax, R1, R2;
    TrivTVStruct *TVs;
    CagdBType
	NewTV = FALSE;

    switch (TV -> GType) {
	case TRIV_TVBEZIER_TYPE:
	    R1 = 0.0;
	    R2 = 1.0;
	    break;
	case TRIV_TVBSPLINE_TYPE:
    	    switch (Dir) {
		case TRIV_CONST_U_DIR:
		    TrivTVDomain(TV, &R1, &R2, &SMin, &SMax, &TMin, &TMax);
		    break;
		case TRIV_CONST_V_DIR:
		    TrivTVDomain(TV, &RMin, &RMax, &R1, &R2, &TMin, &TMax);
		    break;
		case TRIV_CONST_W_DIR:
		    TrivTVDomain(TV, &RMin, &RMax, &SMin, &SMax, &R1, &R2);
		    break;
		default:
		    TRIV_FATAL_ERROR(TRIV_ERR_DIR_NOT_VALID);
		    break;
	    }
	    break;
	default:
	    TRIV_FATAL_ERROR(TRIV_ERR_UNDEF_GEOM);
	    return NULL;
    }

    if (t1 > t2)
	SWAP(CagdRType, t1, t2);

    if (!APX_EQ(t1, R1)) {
	TVs = TrivTVSubdivAtParam(TV, t1, Dir);
	TV = TVs -> Pnext;
	TVs -> Pnext = NULL;
	TrivTVFree(TVs);			   /* Free the first region. */
	NewTV = TRUE;
    }

    if (APX_EQ(t2, R2))
	return NewTV ? TV : TrivTVCopy(TV);
    else {
	TVs = TrivTVSubdivAtParam(TV, t2, Dir);

	if (NewTV)
	    TrivTVFree(TV);

    	TrivTVFree(TVs -> Pnext);		  /* Free the second region. */
    	TVs -> Pnext = NULL;
	return TVs;				/* Returns the first region. */
    }
}

/*****************************************************************************
* DESCRIPTION:                                                               M
* Computes a bounding box for a trivariate freeform function.		     M
*                                                                            *
* PARAMETERS:                                                                M
*   Trivar:   To compute a bounding box for.                                 M
*   BBox:     Where bounding information is to be saved.                     M
*                                                                            *
* RETURN VALUE:                                                              M
*   void                                                                     M
*                                                                            *
* KEYWORDS:                                                                  M
*   TrivTVBBox, bbox, bounding box                                           M
*****************************************************************************/
void TrivTVBBox(TrivTVStruct *Trivar, CagdBBoxStruct *BBox)
{
    TrivTVStruct
	*E3TV = TrivCoerceTVTo(Trivar, CAGD_PT_E3_TYPE);
    int Length = E3TV -> ULength * E3TV -> VLength * E3TV -> WLength;
    CagdRType
	**Points = E3TV -> Points;

    CagdPointsBBox(Points, Length, BBox);

    TrivTVFree(E3TV);
}

/*****************************************************************************
* DESCRIPTION:                                                               M
* Computes a bounding box for a list of trivariate freeform function.        M
*                                                                            *
* PARAMETERS:                                                                M
*   Trivars:  To compute a bounding box for.                                 M
*   BBox:     Where bounding information is to be saved.                     M
*                                                                            *
* RETURN VALUE:                                                              M
*   void                                                                     M
*                                                                            *
* KEYWORDS:                                                                  M
*   TrivTVListBBox, bbox, bounding box                                       M
*****************************************************************************/
void TrivTVListBBox(TrivTVStruct *TVs, CagdBBoxStruct *BBox)
{
    CAGD_RESET_BBOX(BBox);

    for ( ; TVs != NULL; TVs = TVs -> Pnext) {
	CagdBBoxStruct TmpBBox;

	TrivTVBBox(TVs, &TmpBBox);
	CagdMergeBBox(BBox, &TmpBBox);
    }
}
