dnl This is m4 source.
dnl Process using m4 to produce 'C' language file.
dnl
undefine(`begin')dnl
undefine(`index')dnl
undefine(`len')dnl
dnl
dnl If you see this line, you can ignore the next one.
/* Do not edit this file. It is produced from the corresponding .m4 source */
dnl
/*
 *	Copyright 1996, University Corporation for Atmospheric Research
 *      See netcdf/COPYRIGHT file for copying and redistribution conditions.
 */
/* $Id: putget.m4,v 2.62 2006/06/20 15:05:14 ed Exp $ */

#include "nc.h"
#include <string.h>
#include <stdlib.h>
#include <assert.h>
#include "ncx.h"
#include "fbits.h"
#include "onstack.h"
#ifdef LOCKNUMREC
#  include <mpp/shmem.h>	/* for SGI/Cray SHMEM routines */
#  ifdef LN_TEST
#    include <stdio.h>
#  endif
#endif

#undef MIN  /* system may define MIN somewhere and complain */
#define MIN(mm,nn) (((mm) < (nn)) ? (mm) : (nn))

/* #define ODEBUG 1 */

#if ODEBUG
#include <stdio.h>
/*
 * Print the values of an array of size_t
 */
void
arrayp(const char *label, size_t count, const size_t *array)
{
	(void) fprintf(stderr, "%s", label);
	(void) fputc('\t',stderr);	
	for(; count > 0; count--, array++)
		(void) fprintf(stderr," %lu", (unsigned long)*array);
	(void) fputc('\n',stderr);	
}
#endif /* ODEBUG */

/*
 *  This is how much space is required by the user, as in
 *
 *   vals = malloc(nel * nctypelen(var.type));
 *   ncvarget(cdfid, varid, cor, edg, vals);
 */
int
nctypelen(nc_type type) 
{
	switch(type){
	case NC_BYTE :
	case NC_CHAR :
		return((int)sizeof(char));
	case NC_SHORT :
		return(int)(sizeof(short));
	case NC_INT :
		return((int)sizeof(int));
	case NC_FLOAT :
		return((int)sizeof(float));
	case NC_DOUBLE : 
		return((int)sizeof(double));
	}

	return -1;
}


/* Begin fill */
/*
 * This is tunable parameter.
 * It essentially controls the tradeoff between the number of times
 * memcpy() gets called to copy the external data to fill 
 * a large buffer vs the number of times its called to
 * prepare the external data.
 */
#if	_SX
/* NEC SX specific optimization */
#define	NFILL	2048
#else
#define	NFILL	16
#endif


dnl
dnl NCFILL(Type, Xtype, XSize, Fill)
dnl
define(`NCFILL',dnl
`dnl
static int
NC_fill_$2(
	void **xpp,
	size_t nelems)	/* how many */
{
	$1 fillp[NFILL * sizeof(double)/$3];

	assert(nelems <= sizeof(fillp)/sizeof(fillp[0]));

	{
		$1 *vp = fillp;	/* lower bound of area to be filled */
		const $1 *const end = vp + nelems;
		while(vp < end)
		{
			*vp++ = $4;
		}
	}
	return ncx_putn_$2_$1(xpp, nelems, fillp);
}
')dnl

/*
 * Next 6 type specific functions
 * Fill a some memory with the default special value.
 * Formerly
NC_arrayfill()
 */
NCFILL(schar, schar, X_SIZEOF_CHAR, NC_FILL_BYTE)
NCFILL(char, char, X_SIZEOF_CHAR, NC_FILL_CHAR)
NCFILL(short, short, X_SIZEOF_SHORT, NC_FILL_SHORT)

#if (SIZEOF_INT >= X_SIZEOF_INT)
NCFILL(int, int, X_SIZEOF_INT, NC_FILL_INT)
#elif SIZEOF_LONG == X_SIZEOF_INT
NCFILL(long, int, X_SIZEOF_INT, NC_FILL_INT)
#else
#error "NC_fill_int implementation"
#endif

NCFILL(float, float, X_SIZEOF_FLOAT, NC_FILL_FLOAT)
NCFILL(double, double, X_SIZEOF_DOUBLE, NC_FILL_DOUBLE)




/* 
 * Fill the external space for variable 'varp' values at 'recno' with
 * the appropriate value. If 'varp' is not a record variable, fill the
 * whole thing.  For the special case when 'varp' is the only record
 * variable and it is of type byte, char, or short, varsize should be
 * ncp->recsize, otherwise it should be varp->len.
 * Formerly
xdr_NC_fill()
 */
int
fill_NC_var(NC *ncp, const NC_var *varp, size_t varsize, size_t recno)
{
	char xfillp[NFILL * X_SIZEOF_DOUBLE];
	const size_t step = varp->xsz;
	const size_t nelems = sizeof(xfillp)/step;
	const size_t xsz = varp->xsz * nelems;
	NC_attr **attrpp = NULL;
	off_t offset;
	size_t remaining = varsize;

	void *xp;
	int status = NC_NOERR;

	/*
	 * Set up fill value
	 */
	attrpp = NC_findattr(&varp->attrs, _FillValue);
	if( attrpp != NULL )
	{
		/* User defined fill value */
		if( (*attrpp)->type != varp->type || (*attrpp)->nelems != 1 )
		{
			return NC_EBADTYPE;
		}
		else
		{
			/* Use the user defined value */
			char *cp = xfillp;
			const char *const end = &xfillp[sizeof(xfillp)];

			assert(step <= (*attrpp)->xsz);

			for( /*NADA*/; cp < end; cp += step)
			{
				(void) memcpy(cp, (*attrpp)->xvalue, step);
			}
		}
	}
	else
	{
		/* use the default */
		
		assert(xsz % X_ALIGN == 0);
		assert(xsz <= sizeof(xfillp));
	
		xp = xfillp;
	
		switch(varp->type){
		case NC_BYTE :
			status = NC_fill_schar(&xp, nelems);
			break;
		case NC_CHAR :
			status = NC_fill_char(&xp, nelems);
			break;
		case NC_SHORT :
			status = NC_fill_short(&xp, nelems);
			break;
		case NC_INT :
			status = NC_fill_int(&xp, nelems);
			break;
		case NC_FLOAT :
			status = NC_fill_float(&xp, nelems);
			break;
		case NC_DOUBLE : 
			status = NC_fill_double(&xp, nelems);
			break;
		default :
			assert("fill_NC_var invalid type" == 0);
			status = NC_EBADTYPE;
			break;
		}
		if(status != NC_NOERR)
			return status;
	
		assert(xp == xfillp + xsz);
	}

	/*
	 * copyout:
	 * xfillp now contains 'nelems' elements of the fill value
	 * in external representation.
	 */

	/*
	 * Copy it out.
	 */

	offset = varp->begin;
	if(IS_RECVAR(varp))
	{
		offset += (off_t)ncp->recsize * recno;
	}

	assert(remaining > 0);
	for(;;)
	{
		const size_t chunksz = MIN(remaining, ncp->chunk);
		size_t ii;

		status = ncp->nciop->get(ncp->nciop, offset, chunksz,
				 RGN_WRITE, &xp);	
		if(status != NC_NOERR)
		{
			return status;
		}

		/*
		 * fill the chunksz buffer in units  of xsz
		 */
		for(ii = 0; ii < chunksz/xsz; ii++)
		{
			(void) memcpy(xp, xfillp, xsz);
			xp = (char *)xp + xsz;
		}
		/*
		 * Deal with any remainder
		 */
		{
			const size_t rem = chunksz % xsz;
			if(rem != 0)
			{
				(void) memcpy(xp, xfillp, rem);
				/* xp = (char *)xp + xsz; */
			}

		}

		status = ncp->nciop->rel(ncp->nciop, offset, RGN_MODIFIED);

		if(status != NC_NOERR)
		{
			break;
		}

		remaining -= chunksz;
		if(remaining == 0)
			break;	/* normal loop exit */
		offset += chunksz;

	}

	return status;
}
/* End fill */


/*
 * Add a record containing the fill values.
 */
static int
NCfillrecord(NC *ncp, const NC_var *const *varpp, size_t recno)
{
	size_t ii = 0;
	for(; ii < ncp->vars.nelems; ii++, varpp++)
	{
		if( !IS_RECVAR(*varpp) )
		{
			continue;	/* skip non-record variables */
		}
		{
		const int status = fill_NC_var(ncp, *varpp, (*varpp)->len, recno);
		if(status != NC_NOERR)
			return status;
		}
	}
	return NC_NOERR;
}


/*
 * Add a record containing the fill values in the special case when
 * there is exactly one record variable, where we don't require each
 * record to be four-byte aligned (no record padding).
 */
static int
NCfillspecialrecord(NC *ncp, const NC_var *varp, size_t recno)
{
    int status;
    assert(IS_RECVAR(varp));
    status = fill_NC_var(ncp, varp, ncp->recsize, recno);
    if(status != NC_NOERR)
	return status;
    return NC_NOERR;
}


/*
 * It is advantageous to
 * #define TOUCH_LAST
 * when using memory mapped io.
 */
#if TOUCH_LAST
/*
 * Grow the file to a size which can contain recno
 */
static int
NCtouchlast(NC *ncp, const NC_var *const *varpp, size_t recno)
{
	int status = NC_NOERR;
	const NC_var *varp = NULL;
	
	{
	size_t ii = 0;
	for(; ii < ncp->vars.nelems; ii++, varpp++)
	{
		if( !IS_RECVAR(*varpp) )
		{
			continue;	/* skip non-record variables */
		}
		varp = *varpp;
	}
	}
	assert(varp != NULL);
	assert( IS_RECVAR(varp) );
	{
		const off_t offset = varp->begin
				+ (off_t)(recno-1) * (off_t)ncp->recsize
				+ (off_t)(varp->len - varp->xsz);
		void *xp;


		status = ncp->nciop->get(ncp->nciop, offset, varp->xsz,
				 RGN_WRITE, &xp);	
		if(status != NC_NOERR)
			return status;
		(void)memset(xp, 0, varp->xsz);
		status = ncp->nciop->rel(ncp->nciop, offset, RGN_MODIFIED);
	}
	return status;
}
#endif /* TOUCH_LAST */


/*
 * Ensure that the netcdf file has 'numrecs' records,
 * add records and fill as neccessary.
 */
static int
NCvnrecs(NC *ncp, size_t numrecs)
{
	int status = NC_NOERR;
#ifdef LOCKNUMREC
	ushmem_t myticket = 0, nowserving = 0;
	ushmem_t numpe = (ushmem_t) _num_pes();

	/* get ticket and wait */
	myticket = shmem_short_finc((shmem_t *) ncp->lock + LOCKNUMREC_LOCK,
		ncp->lock[LOCKNUMREC_BASEPE]);
#ifdef LN_TEST
		fprintf(stderr,"%d of %d : ticket = %hu\n",
			_my_pe(), _num_pes(), myticket);
#endif
	do {
		shmem_short_get((shmem_t *) &nowserving,
			(shmem_t *) ncp->lock + LOCKNUMREC_SERVING, 1,
			ncp->lock[LOCKNUMREC_BASEPE]);
#ifdef LN_TEST
		fprintf(stderr,"%d of %d : serving = %hu\n",
			_my_pe(), _num_pes(), nowserving);
#endif
		/* work-around for non-unique tickets */
		if (nowserving > myticket && nowserving < myticket + numpe ) {
			/* get a new ticket ... you've been bypassed */ 
			/* and handle the unlikely wrap-around effect */
			myticket = shmem_short_finc(
				(shmem_t *) ncp->lock + LOCKNUMREC_LOCK,
				ncp->lock[LOCKNUMREC_BASEPE]);
#ifdef LN_TEST
				fprintf(stderr,"%d of %d : new ticket = %hu\n",
					_my_pe(), _num_pes(), myticket);
#endif
		}
	} while(nowserving != myticket);
	/* now our turn to check & update value */
#endif

	if(numrecs > NC_get_numrecs(ncp))
	{


#if TOUCH_LAST
		status = NCtouchlast(ncp,
			(const NC_var *const*)ncp->vars.value,
			numrecs);
		if(status != NC_NOERR)
			goto common_return;
#endif /* TOUCH_LAST */

		set_NC_ndirty(ncp);

		if(!NC_dofill(ncp))
		{
			/* Simply set the new numrecs value */
			NC_set_numrecs(ncp, numrecs);
		}
		else
		{
		    /* Treat two cases differently: 
		        - exactly one record variable (no padding)
                        - multiple record variables (each record padded 
                          to 4-byte alignment)
		    */
		    NC_var **vpp = (NC_var **)ncp->vars.value;
		    NC_var *const *const end = &vpp[ncp->vars.nelems];
		    NC_var *recvarp = NULL;	/* last record var */
		    int numrecvars = 0;
		    size_t cur_nrecs;
		    
		    /* determine how many record variables */
		    for( /*NADA*/; vpp < end; vpp++) {
			if(IS_RECVAR(*vpp)) {
			    recvarp = *vpp;
			    numrecvars++;
			}
		    }
		    
		    if (numrecvars != 1) { /* usual case */
			/* Fill each record out to numrecs */
			while((cur_nrecs = NC_get_numrecs(ncp)) < numrecs)
			    {
				status = NCfillrecord(ncp,
					(const NC_var *const*)ncp->vars.value,
					cur_nrecs);
				if(status != NC_NOERR)
				{
					break;
				}
				NC_increase_numrecs(ncp, cur_nrecs +1);
			}
			if(status != NC_NOERR)
				goto common_return;
		    } else {	/* special case */
			/* Fill each record out to numrecs */
			while((cur_nrecs = NC_get_numrecs(ncp)) < numrecs)
			    {
				status = NCfillspecialrecord(ncp,
					recvarp,
					cur_nrecs);
				if(status != NC_NOERR)
				{
					break;
				}
				NC_increase_numrecs(ncp, cur_nrecs +1);
			}
			if(status != NC_NOERR)
				goto common_return;
			
		    }
		}

		if(NC_doNsync(ncp))
		{
			status = write_numrecs(ncp);
		}

	}
common_return:
#ifdef LOCKNUMREC
	/* finished with our lock - increment serving number */
	(void) shmem_short_finc((shmem_t *) ncp->lock + LOCKNUMREC_SERVING,
		ncp->lock[LOCKNUMREC_BASEPE]);
#endif
	return status;
}


/* 
 * Check whether 'coord' values are valid for the variable.
 */
static int
NCcoordck(NC *ncp, const NC_var *varp, const size_t *coord)
{
	const size_t *ip;
	size_t *up;

	if(varp->ndims == 0)
		return NC_NOERR;	/* 'scalar' variable */

	if(IS_RECVAR(varp))
	{
		if(*coord > X_INT_MAX)
			return NC_EINVALCOORDS; /* sanity check */
		if(NC_readonly(ncp) && *coord >= NC_get_numrecs(ncp))
		{
			if(!NC_doNsync(ncp))
				return NC_EINVALCOORDS;
			/* else */
			{
				/* Update from disk and check again */
				const int status = read_numrecs(ncp);
				if(status != NC_NOERR)
					return status;
				if(*coord >= NC_get_numrecs(ncp))
					return NC_EINVALCOORDS;
			}
		}
		ip = coord + 1;
		up = varp->shape + 1;
	}
	else
	{
		ip = coord;
		up = varp->shape;
	}
	
#ifdef CDEBUG
fprintf(stderr,"	NCcoordck: coord %ld, count %d, ip %ld\n",
		coord, varp->ndims, ip );
#endif /* CDEBUG */

	for(; ip < coord + varp->ndims; ip++, up++)
	{

#ifdef CDEBUG
fprintf(stderr,"	NCcoordck: ip %p, *ip %ld, up %p, *up %lu\n",
			ip, *ip, up, *up );
#endif /* CDEBUG */

		/* cast needed for braindead systems with signed size_t */
		if((unsigned long) *ip >= (unsigned long) *up )
			return NC_EINVALCOORDS;
	}

	return NC_NOERR;
}


/* 
 * Check whether 'edges' are valid for the variable and 'start'
 */
/*ARGSUSED*/
static int
NCedgeck(const NC *ncp, const NC_var *varp,
	 const size_t *start, const size_t *edges)
{
	const size_t *const end = start + varp->ndims;
	const size_t *shp = varp->shape;

	if(varp->ndims == 0)
		return NC_NOERR;	/* 'scalar' variable */

	if(IS_RECVAR(varp))
	{
		start++;
		edges++;
		shp++;
	}

	for(; start < end; start++, edges++, shp++)
	{
		/* cast needed for braindead systems with signed size_t */
		if((unsigned long) *edges > *shp ||
			(unsigned long) *start + (unsigned long) *edges > *shp)
		{
			return(NC_EEDGE);
		}
	}
	return NC_NOERR;
}


/* 
 * Translate the (variable, coord) pair into a seek index
 */
static off_t
NC_varoffset(const NC *ncp, const NC_var *varp, const size_t *coord)
{
	if(varp->ndims == 0) /* 'scalar' variable */
		return varp->begin;

	if(varp->ndims == 1)
	{
		if(IS_RECVAR(varp))
			return varp->begin +
				 (off_t)(*coord) * (off_t)ncp->recsize;
		/* else */
		return varp->begin + (off_t)(*coord) * (off_t)varp->xsz;
	}
	/* else */
	{
		off_t lcoord = (off_t)coord[varp->ndims -1];

		size_t *up = varp->dsizes +1;
		const size_t *ip = coord;
		const size_t *const end = varp->dsizes + varp->ndims;
		
		if(IS_RECVAR(varp))
			up++, ip++;

		for(; up < end; up++, ip++)
			lcoord += *up * *ip;

		lcoord *= varp->xsz;
		
		if(IS_RECVAR(varp))
			lcoord += (off_t)(*coord) * ncp->recsize;
		
		lcoord += varp->begin;
		return lcoord;
	}
}


dnl
dnl Output 'nelems' items of contiguous data of type "Type"
dnl for variable 'varp' at 'start'.
dnl "Xtype" had better match 'varp->type'.
dnl---
dnl
dnl PUTNCVX(Xtype, Type)
dnl
define(`PUTNCVX',dnl
`dnl
static int
putNCvx_$1_$2(NC *ncp, const NC_var *varp,
		 const size_t *start, size_t nelems, const $2 *value)
{
	off_t offset = NC_varoffset(ncp, varp, start);
	size_t remaining = varp->xsz * nelems;
	int status = NC_NOERR;
	void *xp;

	if(nelems == 0)
		return NC_NOERR;

	assert(value != NULL);

	for(;;)
	{
		size_t extent = MIN(remaining, ncp->chunk);
		size_t nput = ncx_howmany(varp->type, extent);

		int lstatus = ncp->nciop->get(ncp->nciop, offset, extent,
				 RGN_WRITE, &xp);	
		if(lstatus != NC_NOERR)
			return lstatus;
		
		lstatus = ncx_putn_$1_$2(&xp, nput, value);
		if(lstatus != NC_NOERR && status == NC_NOERR)
		{
			/* not fatal to the loop */
			status = lstatus;
		}

		(void) ncp->nciop->rel(ncp->nciop, offset,
				 RGN_MODIFIED);	

		remaining -= extent;
		if(remaining == 0)
			break; /* normal loop exit */
		offset += extent;
		value += nput;

	}

	return status;
}
')dnl

PUTNCVX(char, char)

PUTNCVX(schar, schar)
PUTNCVX(schar, uchar)
PUTNCVX(schar, short)
PUTNCVX(schar, int)
PUTNCVX(schar, long)
PUTNCVX(schar, float)
PUTNCVX(schar, double)

PUTNCVX(short, schar)
PUTNCVX(short, uchar)
PUTNCVX(short, short)
PUTNCVX(short, int)
PUTNCVX(short, long)
PUTNCVX(short, float)
PUTNCVX(short, double)

PUTNCVX(int, schar)
PUTNCVX(int, uchar)
PUTNCVX(int, short)
PUTNCVX(int, int)
PUTNCVX(int, long)
PUTNCVX(int, float)
PUTNCVX(int, double)

PUTNCVX(float, schar)
PUTNCVX(float, uchar)
PUTNCVX(float, short)
PUTNCVX(float, int)
PUTNCVX(float, long)
PUTNCVX(float, float)
PUTNCVX(float, double)

PUTNCVX(double, schar)
PUTNCVX(double, uchar)
PUTNCVX(double, short)
PUTNCVX(double, int)
PUTNCVX(double, long)
PUTNCVX(double, float)
PUTNCVX(double, double)


dnl
dnl PUTNCV(Type)
dnl
define(`PUTNCV',dnl
`dnl
static int
putNCv_$1(NC *ncp, const NC_var *varp,
		 const size_t *start, size_t nelems, const $1 *value)
{
	switch(varp->type){
	case NC_CHAR:
		return NC_ECHAR;
	case NC_BYTE:
		return putNCvx_schar_$1(ncp, varp, start, nelems,
			value);
	case NC_SHORT:
		return putNCvx_short_$1(ncp, varp, start, nelems,
			value);
	case NC_INT:
		return putNCvx_int_$1(ncp, varp, start, nelems,
			value);
	case NC_FLOAT:
		return putNCvx_float_$1(ncp, varp, start, nelems,
			value);
	case NC_DOUBLE: 
		return putNCvx_double_$1(ncp, varp, start, nelems,
			value);
	}
	return NC_EBADTYPE;
}
')dnl

static int
putNCv_text(NC *ncp, const NC_var *varp,
		 const size_t *start, size_t nelems, const char *value)
{
	if(varp->type != NC_CHAR)
		return NC_ECHAR;
	return putNCvx_char_char(ncp, varp, start, nelems, value);
}

PUTNCV(schar)
PUTNCV(uchar)
PUTNCV(short)
PUTNCV(int)
PUTNCV(long)
PUTNCV(float)
PUTNCV(double)


dnl
dnl GETNCVX(XType, Type)
dnl
define(`GETNCVX',dnl
`dnl
static int
getNCvx_$1_$2(const NC *ncp, const NC_var *varp,
		 const size_t *start, size_t nelems, $2 *value)
{
	off_t offset = NC_varoffset(ncp, varp, start);
	size_t remaining = varp->xsz * nelems;
	int status = NC_NOERR;
	const void *xp;

	if(nelems == 0)
		return NC_NOERR;

	assert(value != NULL);

	for(;;)
	{
		size_t extent = MIN(remaining, ncp->chunk);
		size_t nget = ncx_howmany(varp->type, extent);

		int lstatus = ncp->nciop->get(ncp->nciop, offset, extent,
				 0, (void **)&xp);	/* cast away const */
		if(lstatus != NC_NOERR)
			return lstatus;
		
		lstatus = ncx_getn_$1_$2(&xp, nget, value);
		if(lstatus != NC_NOERR && status == NC_NOERR)
			status = lstatus;

		(void) ncp->nciop->rel(ncp->nciop, offset, 0);	

		remaining -= extent;
		if(remaining == 0)
			break; /* normal loop exit */
		offset += extent;
		value += nget;
	}

	return status;
}
')dnl

GETNCVX(char, char)

GETNCVX(schar, schar)
GETNCVX(schar, uchar)
GETNCVX(schar, short)
GETNCVX(schar, int)
GETNCVX(schar, long)
GETNCVX(schar, float)
GETNCVX(schar, double)

GETNCVX(short, schar)
GETNCVX(short, uchar)
GETNCVX(short, short)
GETNCVX(short, int)
GETNCVX(short, long)
GETNCVX(short, float)
GETNCVX(short, double)

GETNCVX(int, schar)
GETNCVX(int, uchar)
GETNCVX(int, short)
GETNCVX(int, int)
GETNCVX(int, long)
GETNCVX(int, float)
GETNCVX(int, double)

GETNCVX(float, schar)
GETNCVX(float, uchar)
GETNCVX(float, short)
GETNCVX(float, int)
GETNCVX(float, long)
GETNCVX(float, float)
GETNCVX(float, double)

GETNCVX(double, schar)
GETNCVX(double, uchar)
GETNCVX(double, short)
GETNCVX(double, int)
GETNCVX(double, long)
GETNCVX(double, float)
GETNCVX(double, double)


dnl
dnl GETNCV(Type)
dnl
define(`GETNCV',dnl
`dnl
static int
getNCv_$1(const NC *ncp, const NC_var *varp,
		 const size_t *start, size_t nelems, $1 *value)
{
	switch(varp->type){
	case NC_CHAR:
		return NC_ECHAR;
	case NC_BYTE:
		return getNCvx_schar_$1(ncp, varp, start, nelems,
			value);
	case NC_SHORT:
		return getNCvx_short_$1(ncp, varp, start, nelems,
			value);
	case NC_INT:
		return getNCvx_int_$1(ncp, varp, start, nelems,
			value);
	case NC_FLOAT:
		return getNCvx_float_$1(ncp, varp, start, nelems,
			value);
	case NC_DOUBLE: 
		return getNCvx_double_$1(ncp, varp, start, nelems,
			value);
	}
	return NC_EBADTYPE;
}
')dnl

GETNCV(schar)
GETNCV(uchar)
GETNCV(short)
GETNCV(int)
GETNCV(long)
GETNCV(float)
GETNCV(double)


static int
getNCv_text(const NC *ncp, const NC_var *varp,
		 const size_t *start, size_t nelems, char *value)
{
	if(varp->type != NC_CHAR)
		return NC_ECHAR;
	return getNCvx_char_char(ncp, varp, start, nelems, value);
}


/*
 * Copy 'nbytes' contiguous external values
 * from ('inncp', invp', inncoord')
 * to   ('outncp', 'outvp', 'outcoord')
 * 'inncp' shouldn't be the same as 'outncp'.
 * Used only by ncvarcopy()
 */
static int
NCxvarcpy(NC *inncp, NC_var *invp, size_t *incoord,
	NC *outncp, NC_var *outvp, size_t *outcoord, size_t nbytes)
{
	int status;
	off_t inoffset = NC_varoffset(inncp, invp, incoord);
	off_t outoffset = NC_varoffset(outncp, outvp, outcoord);
	void *inxp;
	void *outxp;
	const size_t chunk = MIN(inncp->chunk, outncp->chunk);

	do {
		const size_t extent = MIN(nbytes, chunk);

		status = inncp->nciop->get(inncp->nciop, inoffset, extent,
				 0, &inxp);	
		if(status != NC_NOERR)
			return status;

		status = outncp->nciop->get(outncp->nciop, outoffset, extent,
				 RGN_WRITE, &outxp);	
		if(status != NC_NOERR)
		{
			(void) inncp->nciop->rel(inncp->nciop, inoffset, 0);	
			break;
		}

		(void) memcpy(outxp, inxp, extent);

		status = outncp->nciop->rel(outncp->nciop, outoffset,
			 RGN_MODIFIED);
		(void) inncp->nciop->rel(inncp->nciop, inoffset, 0);	

		nbytes -= extent;
		if(nbytes == 0)
			break; /* normal loop exit */
		inoffset += extent;
		outoffset += extent;
		
	} while (status == NC_NOERR);

	return status;
}


/*
 *  For ncvar{put,get},
 *  find the largest contiguous block from within 'edges'.
 *  returns the index to the left of this (which may be -1).
 *  Compute the number of contiguous elements and return
 *  that in *iocountp.
 *  The presence of "record" variables makes this routine
 *  overly subtle.
 */
static int
NCiocount(const NC *const ncp, const NC_var *const varp,
	const size_t *const edges,
	size_t *const iocountp)
{
	const size_t *edp0 = edges;
	const size_t *edp = edges + varp->ndims;
	const size_t *shp = varp->shape + varp->ndims;

	if(IS_RECVAR(varp))
	{
		if(varp->ndims == 1 && ncp->recsize <= varp->len)
		{
			/* one dimensional && the only 'record' variable */
			*iocountp = *edges;
			return(0);
		}
		/* else */
		edp0++;
	}

	assert(edges != NULL);

	/* find max contiguous */
	while(edp > edp0)
	{
		shp--; edp--;
		if(*edp < *shp )
		{
			const size_t *zedp = edp;
			while(zedp >= edp0)
			{
				if(*zedp == 0)
				{
					*iocountp = 0;
					goto done;
				}
				/* Tip of the hat to segmented architectures */
				if(zedp == edp0)
					break;
				zedp--;
			}
			break;
		}
		assert(*edp == *shp);
	}

	/*
	 * edp, shp reference rightmost index s.t. *(edp +1) == *(shp +1)
	 *
	 * Or there is only one dimension.
	 * If there is only one dimension and it is 'non record' dimension,
	 * 	edp is &edges[0] and we will return -1.
	 * If there is only one dimension and and it is a "record dimension",
	 *	edp is &edges[1] (out of bounds) and we will return 0;
	 */
	assert(shp >= varp->shape + varp->ndims -1 
		|| *(edp +1) == *(shp +1));

	/* now accumulate max count for a single io operation */
	for(*iocountp = 1, edp0 = edp;
		 	edp0 < edges + varp->ndims;
			edp0++)
	{
		*iocountp *= *edp0;
	}

done:
	return((int)(edp - edges) - 1);
}


/*
 * Set the elements of the array 'upp' to
 * the sum of the corresponding elements of
 * 'stp' and 'edp'. 'end' should be &stp[nelems].
 */
static void
set_upper(size_t *upp, /* modified on return */
	const size_t *stp,
	const size_t *edp,
	const size_t *const end)
{
	while(upp < end) {
		*upp++ = *stp++ + *edp++;
	}
}


/*
 * The infamous and oft-discussed odometer code.
 *
 * 'start[]' is the starting coordinate.
 * 'upper[]' is the upper bound s.t. start[ii] < upper[ii].
 * 'coord[]' is the register, the current coordinate value.
 * For some ii,
 * upp == &upper[ii]
 * cdp == &coord[ii]
 * 
 * Running this routine increments *cdp.
 *
 * If after the increment, *cdp is equal to *upp
 * (and cdp is not the leftmost dimension),
 * *cdp is "zeroed" to the starting value and
 * we need to "carry", eg, increment one place to
 * the left.
 * 
 * TODO: Some architectures hate recursion?
 * 	Reimplement non-recursively.
 */
static void
odo1(const size_t *const start, const size_t *const upper,
	size_t *const coord, /* modified on return */
	const size_t *upp,
	size_t *cdp)
{
	assert(coord <= cdp && cdp <= coord + NC_MAX_VAR_DIMS);
	assert(upper <= upp && upp <= upper + NC_MAX_VAR_DIMS);
	assert(upp - upper == cdp - coord);
	
	assert(*cdp <= *upp);

	(*cdp)++;
	if(cdp != coord && *cdp >= *upp)
	{
		*cdp = start[cdp - coord];
		odo1(start, upper, coord, upp -1, cdp -1);
	}
}
#ifdef _CRAYC
#pragma _CRI noinline odo1
#endif


dnl
dnl NCTEXTCOND(Abbrv)
dnl This is used inside the NC{PUT,GET} macros below
dnl
define(`NCTEXTCOND',dnl
`dnl
ifelse($1, text,dnl
`dnl
	if(varp->type != NC_CHAR)
		return NC_ECHAR;
',dnl
`dnl
	if(varp->type == NC_CHAR)
		return NC_ECHAR;
')dnl
')dnl


/* Public */

dnl
dnl NCPUTVAR1(Abbrev, Type)
dnl
define(`NCPUTVAR1',dnl
`dnl
int
nc_put_var1_$1(int ncid, int varid, const size_t *coord,
	const $2 *value)
{
	int status;
	NC *ncp;
	const NC_var *varp;

	status = NC_check_id(ncid, &ncp); 
	if(status != NC_NOERR)
		return status;

	if(NC_readonly(ncp))
		return NC_EPERM;

	if(NC_indef(ncp))
		return NC_EINDEFINE;

	varp = NC_lookupvar(ncp, varid);
	if(varp == NULL)
		return NC_ENOTVAR; /* TODO: lost NC_EGLOBAL */

NCTEXTCOND($1)
	status = NCcoordck(ncp, varp, coord);
	if(status != NC_NOERR)
		return status;

	if(IS_RECVAR(varp))
	{
		status = NCvnrecs(ncp, *coord +1);
		if(status != NC_NOERR)
			return status;
	}

	return putNCv_$1(ncp, varp, coord, 1, value);
}
')dnl

NCPUTVAR1(text, char)

NCPUTVAR1(uchar, uchar)
NCPUTVAR1(schar, schar)
NCPUTVAR1(short, short)
NCPUTVAR1(int, int)
NCPUTVAR1(long, long)
NCPUTVAR1(float, float)
NCPUTVAR1(double, double)

dnl
dnl NCGETVAR1(Abbrv, Type)
dnl
define(`NCGETVAR1',dnl
`dnl
int
nc_get_var1_$1(int ncid, int varid, const size_t *coord, $2 *value)
{
	int status;
	NC *ncp;
	const NC_var *varp;

	status = NC_check_id(ncid, &ncp); 
	if(status != NC_NOERR)
		return status;

	if(NC_indef(ncp))
		return NC_EINDEFINE;

	varp = NC_lookupvar(ncp, varid);
	if(varp == NULL)
		return NC_ENOTVAR; /* TODO: lost NC_EGLOBAL */

NCTEXTCOND($1)
	status = NCcoordck(ncp, varp, coord);
	if(status != NC_NOERR)
		return status;

	return getNCv_$1(ncp, varp, coord, 1, value);
}
')dnl

NCGETVAR1(text, char)

NCGETVAR1(uchar, uchar)
NCGETVAR1(schar, schar)
NCGETVAR1(short, short)
NCGETVAR1(int, int)
NCGETVAR1(long, long)
NCGETVAR1(float, float)
NCGETVAR1(double, double)

dnl
dnl NCPUTVARA(Abbrv, Type)
dnl
define(`NCPUTVARA',dnl
`dnl
int
nc_put_vara_$1(int ncid, int varid,
	 const size_t *start, const size_t *edges, const $2 *value)
{
	int status = NC_NOERR;
	NC *ncp;
	const NC_var *varp;
	int ii;
	size_t iocount;

	status = NC_check_id(ncid, &ncp); 
	if(status != NC_NOERR)
		return status;

	if(NC_readonly(ncp))
		return NC_EPERM;

	if(NC_indef(ncp))
		return NC_EINDEFINE;

	varp = NC_lookupvar(ncp, varid);
	if(varp == NULL)
		return NC_ENOTVAR; /* TODO: lost NC_EGLOBAL */

NCTEXTCOND($1)
	status = NCcoordck(ncp, varp, start);
	if(status != NC_NOERR)
		return status;
	status = NCedgeck(ncp, varp, start, edges);
	if(status != NC_NOERR)
		return status;

	if(varp->ndims == 0) /* scalar variable */
	{
		return( putNCv_$1(ncp, varp, start, 1, value) );
	}

	if(IS_RECVAR(varp))
	{
		status = NCvnrecs(ncp, *start + *edges);
		if(status != NC_NOERR)
			return status;

		if(varp->ndims == 1
			&& ncp->recsize <= varp->len)
		{
			/* one dimensional && the only record variable  */
			return( putNCv_$1(ncp, varp, start, *edges, value) );
		}
	}

	/*
	 * find max contiguous
	 *   and accumulate max count for a single io operation
	 */
	ii = NCiocount(ncp, varp, edges, &iocount);

	if(ii == -1)
	{
		return( putNCv_$1(ncp, varp, start, iocount, value) );
	}

	assert(ii >= 0);


	{ /* inline */
	ALLOC_ONSTACK(coord, size_t, varp->ndims);
	ALLOC_ONSTACK(upper, size_t, varp->ndims);
	const size_t index = ii;

	/* copy in starting indices */
	(void) memcpy(coord, start, varp->ndims * sizeof(size_t));

	/* set up in maximum indices */
	set_upper(upper, start, edges, &upper[varp->ndims]);

	/* ripple counter */
	while(*coord < *upper)
	{
		const int lstatus = putNCv_$1(ncp, varp, coord, iocount,
				 value);
		if(lstatus != NC_NOERR)
		{
			if(lstatus != NC_ERANGE)
			{
				status = lstatus;
				/* fatal for the loop */
				break;
			}
			/* else NC_ERANGE, not fatal for the loop */
			if(status == NC_NOERR)
				status = lstatus;
		}
		value += iocount;
		odo1(start, upper, coord, &upper[index], &coord[index]);
	}

	FREE_ONSTACK(upper);
	FREE_ONSTACK(coord);
	} /* end inline */

	return status;
}
')dnl

NCPUTVARA(text, char)

NCPUTVARA(uchar, uchar)
NCPUTVARA(schar, schar)
NCPUTVARA(short, short)
NCPUTVARA(int, int)
NCPUTVARA(long, long)
NCPUTVARA(float, float)
NCPUTVARA(double, double)


dnl
dnl NCGETVARA(Abbrv, Type)
dnl
define(`NCGETVARA',dnl
`dnl
int
nc_get_vara_$1(int ncid, int varid,
	 const size_t *start, const size_t *edges, $2 *value)
{
	int status = NC_NOERR;
	NC *ncp;
	const NC_var *varp;
	int ii;
	size_t iocount;

	status = NC_check_id(ncid, &ncp); 
	if(status != NC_NOERR)
		return status;

	if(NC_indef(ncp))
		return NC_EINDEFINE;

	varp = NC_lookupvar(ncp, varid);
	if(varp == NULL)
		return NC_ENOTVAR; /* TODO: lost NC_EGLOBAL */

NCTEXTCOND($1)
	status = NCcoordck(ncp, varp, start);
	if(status != NC_NOERR)
		return status;
	status = NCedgeck(ncp, varp, start, edges);
	if(status != NC_NOERR)
		return status;

	if(varp->ndims == 0) /* scalar variable */
	{
		return( getNCv_$1(ncp, varp, start, 1, value) );
	}

	if(IS_RECVAR(varp))
	{
		if(*start + *edges > NC_get_numrecs(ncp))
			return NC_EEDGE;
		if(varp->ndims == 1 && ncp->recsize <= varp->len)
		{
			/* one dimensional && the only record variable  */
			return( getNCv_$1(ncp, varp, start, *edges, value) );
		}
	}

	/*
	 * find max contiguous
	 *   and accumulate max count for a single io operation
	 */
	ii = NCiocount(ncp, varp, edges, &iocount);

	if(ii == -1)
	{
		return( getNCv_$1(ncp, varp, start, iocount, value) );
	}

	assert(ii >= 0);


	{ /* inline */
	ALLOC_ONSTACK(coord, size_t, varp->ndims);
	ALLOC_ONSTACK(upper, size_t, varp->ndims);
	const size_t index = ii;

	/* copy in starting indices */
	(void) memcpy(coord, start, varp->ndims * sizeof(size_t));

	/* set up in maximum indices */
	set_upper(upper, start, edges, &upper[varp->ndims]);

	/* ripple counter */
	while(*coord < *upper)
	{
		const int lstatus = getNCv_$1(ncp, varp, coord, iocount,
				value);
		if(lstatus != NC_NOERR)
		{
			if(lstatus != NC_ERANGE)
			{
				status = lstatus;
				/* fatal for the loop */
				break;
			}
			/* else NC_ERANGE, not fatal for the loop */
			if(status == NC_NOERR)
				status = lstatus;
		}
		value += iocount;
		odo1(start, upper, coord, &upper[index], &coord[index]);
	}

	FREE_ONSTACK(upper);
	FREE_ONSTACK(coord);
	} /* end inline */

	return status;
}
')dnl

NCGETVARA(text, char)

NCGETVARA(uchar, uchar)
NCGETVARA(schar, schar)
NCGETVARA(short, short)
NCGETVARA(int, int)
NCGETVARA(long, long)
NCGETVARA(float, float)
NCGETVARA(double, double)


#if defined(__cplusplus)
/* C++ consts default to internal linkage and must be initialized */
const size_t coord_zero[NC_MAX_VAR_DIMS] = {0};
#else
static const size_t coord_zero[NC_MAX_VAR_DIMS];
#endif

dnl
dnl NCPUTVAR(Abbrev, Type)
dnl
define(`NCPUTVAR',dnl
`dnl
int
nc_put_var_$1(int ncid, int varid, const $2 *value)
{
	int status = NC_NOERR;
	NC *ncp;
	const NC_var *varp;

	status = NC_check_id(ncid, &ncp); 
	if(status != NC_NOERR)
		return status;

	if(NC_readonly(ncp))
		return NC_EPERM;

	if(NC_indef(ncp))
		return NC_EINDEFINE;

	varp = NC_lookupvar(ncp, varid);
	if(varp == NULL)
		return NC_ENOTVAR; /* TODO: lost NC_EGLOBAL */

NCTEXTCOND($1)
	if(varp->ndims == 0) /* scalar variable */
	{
		const size_t zed = 0;
		return( putNCv_$1(ncp, varp, &zed, 1, value) );
	}

	if(!IS_RECVAR(varp))
	{
		return(putNCv_$1(ncp, varp, coord_zero, *varp->dsizes, value));
	}
	/* else */

	if(varp->ndims == 1
			&& ncp->recsize <= varp->len)
	{
		/* one dimensional && the only record variable  */
		return(putNCv_$1(ncp, varp, coord_zero, NC_get_numrecs(ncp),
			value));
	}
	/* else */

	{
	ALLOC_ONSTACK(coord, size_t, varp->ndims);
	size_t elemsPerRec = 1;
	const size_t nrecs = NC_get_numrecs(ncp);
	(void) memset(coord, 0, varp->ndims * sizeof(size_t));
	/* TODO: fix dsizes to avoid this nonsense */
	if(varp->ndims > 1)
		elemsPerRec = varp->dsizes[1];
	while(*coord < nrecs)
	{
		const int lstatus = putNCv_$1(ncp, varp, coord, elemsPerRec,
				 value);
		if(lstatus != NC_NOERR)
		{
			if(lstatus != NC_ERANGE)
			{
				status = lstatus;
				/* fatal for the loop */
				break;
			}
			/* else NC_ERANGE, not fatal for the loop */
			if(status == NC_NOERR)
				status = lstatus;
		}
		value += elemsPerRec;
		(*coord)++;
	}
	FREE_ONSTACK(coord);
	} /* elemsPerRec */

	return status;
}
')dnl

NCPUTVAR(text, char)

NCPUTVAR(uchar, uchar)
NCPUTVAR(schar, schar)
NCPUTVAR(short, short)
NCPUTVAR(int, int)
NCPUTVAR(long, long)
NCPUTVAR(float, float)
NCPUTVAR(double, double)


dnl
dnl NCGETVAR(Abbrv, Type)
dnl
define(`NCGETVAR',dnl
`dnl
int
nc_get_var_$1(int ncid, int varid, $2 *value)
{
	int status = NC_NOERR;
	NC *ncp;
	const NC_var *varp;

	status = NC_check_id(ncid, &ncp); 
	if(status != NC_NOERR)
		return status;

	if(NC_indef(ncp))
		return NC_EINDEFINE;

	varp = NC_lookupvar(ncp, varid);
	if(varp == NULL)
		return NC_ENOTVAR; /* TODO: lost NC_EGLOBAL */

	if(varp->ndims == 0) /* scalar variable */
	{
		const size_t zed = 0;
		return( getNCv_$1(ncp, varp, &zed, 1, value) );
	}

NCTEXTCOND($1)

	if(!IS_RECVAR(varp))
	{
		return(getNCv_$1(ncp, varp, coord_zero, *varp->dsizes, value));
	}
	/* else */

	if(varp->ndims == 1
			&& ncp->recsize <= varp->len)
	{
		/* one dimensional && the only record variable  */
		return(getNCv_$1(ncp, varp, coord_zero, NC_get_numrecs(ncp),
			value));
	}
	/* else */

	{
	ALLOC_ONSTACK(coord, size_t, varp->ndims);
	size_t elemsPerRec = 1;
	const size_t nrecs = NC_get_numrecs(ncp);
	(void) memset(coord, 0, varp->ndims * sizeof(size_t));
	/* TODO: fix dsizes to avoid this nonsense */
	if(varp->ndims > 1)
		elemsPerRec = varp->dsizes[1];
	while(*coord < nrecs)
	{
		const int lstatus = getNCv_$1(ncp, varp, coord, elemsPerRec,
				value);
		if(lstatus != NC_NOERR)
		{
			if(lstatus != NC_ERANGE)
			{
				status = lstatus;
				/* fatal for the loop */
				break;
			}
			/* else NC_ERANGE, not fatal for the loop */
			if(status == NC_NOERR)
				status = lstatus;
		}
		value += elemsPerRec;
		(*coord)++;
	}
	FREE_ONSTACK(coord);
	} /* elemsPerRec */

	return status;
}
')dnl

NCGETVAR(text, char)

NCGETVAR(uchar, uchar)
NCGETVAR(schar, schar)
NCGETVAR(short, short)
NCGETVAR(int, int)
NCGETVAR(long, long)
NCGETVAR(float, float)
NCGETVAR(double, double)


/* Begin putgetg.c */

dnl
dnl  NC_VARM_Upper_Body(void)
dnl
define(`NC_VARM_Upper_Body',dnl
`dnl
	int status = ENOERR;
	NC *ncp;
	NC_var *varp;
	int maxidim;	/* maximum dimensional index */

	status = NC_check_id (ncid, &ncp);
	if (status != NC_NOERR)
		return status;

	if (NC_indef (ncp))
	{
		return NC_EINDEFINE;
	}
')dnl
dnl
dnl  NC_VARM_Mid_Body(put_or_get, Abbrv)
dnl
define(`NC_VARM_Mid_Body',dnl
`dnl
	varp = NC_lookupvar (ncp, varid);
	if (varp == NULL)
		return NC_ENOTVAR;

NCTEXTCOND($2)
	maxidim = (int) varp->ndims - 1;

	if (maxidim < 0)
	{
		/*
		 * The variable is a scalar; consequently,
		 * there s only one thing to get and only one place to put it.
		 * (Why was I called?)
		 */
		return $1 (ncp, varp, start, 1, value);
	}
	
	/*
	 * else
	 * The variable is an array.
	 */
	{
		int idim;
		size_t *mystart = NULL;
		size_t *myedges;
		size_t *iocount;	/* count vector */
		size_t *stop;	/* stop indexes */
		size_t *length;	/* edge lengths in bytes */
		ptrdiff_t *mystride;
		ptrdiff_t *mymap;

		/*
		 * Verify stride argument.
		 */
		for (idim = 0; idim <= maxidim; ++idim)
		{
			if (stride != NULL
				&& (stride[idim] == 0
		/* cast needed for braindead systems with signed size_t */
				|| (unsigned long) stride[idim] >= X_INT_MAX))
			{
				return NC_ESTRIDE;
			}
		}

		/* assert(sizeof(ptrdiff_t) >= sizeof(size_t)); */
		mystart = (size_t *)calloc(varp->ndims * 7, sizeof(ptrdiff_t));
		if(mystart == NULL)
			return NC_ENOMEM;
		myedges = mystart + varp->ndims;
		iocount = myedges + varp->ndims;
		stop = iocount + varp->ndims;
		length = stop + varp->ndims;
		mystride = (ptrdiff_t *)(length + varp->ndims);
		mymap = mystride + varp->ndims;

		/*
		 * Initialize I/O parameters.
		 */
		for (idim = maxidim; idim >= 0; --idim)
		{
			mystart[idim] = start != NULL
				? start[idim]
				: 0;

			if (edges[idim] == 0)
			{
				status = NC_NOERR;	/* read/write no data */
				goto done;
			}

			myedges[idim] = edges != NULL
				? edges[idim]
				: idim == 0 && IS_RECVAR (varp)
				? NC_get_numrecs(ncp) - mystart[idim]
				: varp->shape[idim] - mystart[idim];
			mystride[idim] = stride != NULL
				? stride[idim]
				: 1;
			mymap[idim] = map != NULL
				? map[idim]
				: idim == maxidim
				? 1
				: mymap[idim + 1] * (ptrdiff_t) myedges[idim + 1];

			iocount[idim] = 1;
			length[idim] = mymap[idim] * myedges[idim];
			stop[idim] = mystart[idim] + myedges[idim] * mystride[idim];
		}
')dnl
dnl
dnl  NC_VARM_Lower_Body(put_or_get)
dnl
define(`NC_VARM_Lower_Body',dnl
`dnl
		/*
		 * As an optimization, adjust I/O parameters when the fastest 
		 * dimension has unity stride both externally and internally.
		 * In this case, the user could have called a simpler routine
		 * (i.e. ncvar$1()
		 */
		if (mystride[maxidim] == 1
			&& mymap[maxidim] == 1)
		{
			iocount[maxidim] = myedges[maxidim];
			mystride[maxidim] = (ptrdiff_t) myedges[maxidim];
			mymap[maxidim] = (ptrdiff_t) length[maxidim];
		}

		/*
		 * Perform I/O.  Exit when done.
		 */
		for (;;)
		{
			/* TODO: */
			int lstatus = $1 (ncid, varid, mystart, iocount,
						value);
			if (lstatus != NC_NOERR 
				&& (status == NC_NOERR || lstatus != NC_ERANGE))
				status = lstatus;

			/*
			 * The following code permutes through the variable s
			 * external start-index space and it s internal address
			 * space.  At the UPC, this algorithm is commonly
			 * called "odometer code".
			 */
			idim = maxidim;
		carry:
			value += mymap[idim];
			mystart[idim] += mystride[idim];
			if (mystart[idim] == stop[idim])
			{
				mystart[idim] = start[idim];
				value -= length[idim];
				if (--idim < 0)
					break; /* normal return */
				goto carry;
			}
		} /* I/O loop */
	done:
		free(mystart);
	} /* variable is array */
	return status;
')dnl

dnl
dnl NCGETVARS(Abbrv, Type)
dnl
define(`NCGETVARS',dnl
`dnl
int
nc_get_vars_$1 (
	int ncid,
	int varid,
	const size_t * start,
	const size_t * edges,
	const ptrdiff_t * stride,
	$2 *value)
{
	return nc_get_varm_$1 (ncid, varid, start, edges,
			 stride, 0, value);
}
')dnl

NCGETVARS(text, char)

NCGETVARS(uchar, uchar)
NCGETVARS(schar, schar)
NCGETVARS(short, short)
NCGETVARS(int, int)
NCGETVARS(long, long)
NCGETVARS(float, float)
NCGETVARS(double, double)


dnl
dnl NCPUTVARS(Abbrv, Type)
dnl
define(`NCPUTVARS',dnl
`dnl
int
nc_put_vars_$1 (
	int ncid,
	int varid,
	const size_t * start,
	const size_t * edges,
	const ptrdiff_t * stride,
	const $2 *value)
{
	return nc_put_varm_$1 (ncid, varid, start, edges,
			 stride, 0, value);
}
')dnl

NCPUTVARS(text, char)

NCPUTVARS(uchar, uchar)
NCPUTVARS(schar, schar)
NCPUTVARS(short, short)
NCPUTVARS(int, int)
NCPUTVARS(long, long)
NCPUTVARS(float, float)
NCPUTVARS(double, double)


/*
 * Generalized hyperslab input.
 */
dnl
dnl NCGETVARM(Abbrv, Type)
dnl
define(`NCGETVARM',dnl
`dnl
int
nc_get_varm_$1(int ncid, int varid,
	const size_t *start, const size_t *edges,
	const ptrdiff_t *stride,
	const ptrdiff_t *map,
	$2 *value)
{
NC_VARM_Upper_Body()
NC_VARM_Mid_Body(getNCv_$1, $1)
		/*
		 * Check start, edges
		 */
		for (idim = maxidim; idim >= 0; --idim)
		{
			size_t dimlen = 
				idim == 0 && IS_RECVAR (varp)
					? NC_get_numrecs(ncp)
					  : varp->shape[idim];
			if (mystart[idim] >= dimlen)
			{
				status = NC_EINVALCOORDS;
				goto done;
			}

			if (mystart[idim] + myedges[idim] > dimlen)
			{
				status = NC_EEDGE;
				goto done;
			}

		}
NC_VARM_Lower_Body(nc_get_vara_$1)
}
')dnl

NCGETVARM(text, char)

NCGETVARM(uchar, uchar)
NCGETVARM(schar, schar)
NCGETVARM(short, short)
NCGETVARM(int, int)
NCGETVARM(long, long)
NCGETVARM(float, float)
NCGETVARM(double, double)

#ifdef NO_NETCDF_2
extern int
nctypelen(nc_type datatype);
#endif


/*
 * Generalized hyperslab output.
 */
dnl
dnl NCPUTVARM(Abbrv, Type)
dnl
define(`NCPUTVARM',dnl
`dnl
int
nc_put_varm_$1(int ncid, int varid,
	const size_t *start, const size_t *edges,
	const ptrdiff_t *stride, const ptrdiff_t *map,
	const $2 *value)
{
NC_VARM_Upper_Body()
		if (NC_readonly (ncp))
			return NC_EPERM;
NC_VARM_Mid_Body(putNCv_$1, $1)
		/*
		 * Check start, edges
		 */
		for (idim = IS_RECVAR (varp); idim < maxidim; ++idim)
		{
			if (mystart[idim] > varp->shape[idim])
			{
				status = NC_EINVALCOORDS;
				goto done;
			}
			if (mystart[idim] + myedges[idim] > varp->shape[idim])
			{
				status = NC_EEDGE;
				goto done;
			}
		}
NC_VARM_Lower_Body(nc_put_vara_$1)
}
')dnl

NCPUTVARM(text, char)

NCPUTVARM(uchar, uchar)
NCPUTVARM(schar, schar)
NCPUTVARM(short, short)
NCPUTVARM(int, int)
NCPUTVARM(long, long)
NCPUTVARM(float, float)
NCPUTVARM(double, double)


/*
 * Copy the values of a variable from an input netCDF to an output netCDF.
 * Input and output var assummed to have the same shape.
 * return -1 on error.
 */
int
nc_copy_var(int ncid_in, int varid, int ncid_out)
{
	int status = NC_NOERR;
	NC *inncp, *outncp;
	NC_var *invp, *outvp;

	status = NC_check_id(ncid_in, &inncp); 
	if(status != NC_NOERR)
		return status;


	if(NC_indef(inncp))
	{
		return NC_EINDEFINE;
	}

	status = NC_check_id(ncid_out, &outncp); 
	if(status != NC_NOERR)
		return status;

	if(NC_readonly(outncp))
	{
		/* output file isn't writable */
		return NC_EPERM;
	}

	if(NC_indef(outncp))
	{
		return NC_EINDEFINE;
	}

	/* find the variable in the input cdf */
	invp = NC_lookupvar(inncp, varid);
	if(invp == NULL)
	{
		return NC_ENOTVAR;
	}

	/* find the variable in the output cdf */
	if(NC_findvar(&outncp->vars, invp->name->cp, &outvp) == -1)
	{
		return NC_ENOTVAR;
	}

	/* can we even attempt to copy without conversion? */
	if(outvp->type != invp->type)
	{
		return NC_EINVAL;
	}

	if(        (invp->ndims == 0 && outvp->ndims != 0)
		|| (invp->ndims != 0 && outvp->ndims == 0)
		|| (IS_RECVAR(invp) && !IS_RECVAR(outvp))
		|| (!IS_RECVAR(invp) && IS_RECVAR(outvp))
		|| (invp->len != outvp->len)
	)
	{
		return NC_EINVAL;
	}

	/*
	 * Check coordinates
	 */
	{
	ALLOC_ONSTACK(coord, size_t, invp->ndims);
	const size_t nrecs = NC_get_numrecs(inncp);
	(void) memcpy(coord, invp->shape, invp->ndims * sizeof(size_t));
	if(IS_RECVAR(invp))
		*coord = nrecs;
	
	{
	size_t ii = 0;
	for(; ii < invp->ndims; ii++)
		coord[ii] --;
	}
	/* at this point, coord is the largest valid coord of invp */

	if(NCcoordck(outncp, outvp, coord) != NC_NOERR)
	{
		return NC_EINVAL;
	}
	/* else */

	(void) memset(coord, 0, invp->ndims * sizeof(size_t));
	
	if(!IS_RECVAR(invp))
	{
		status = NCxvarcpy(inncp, invp, coord,
				outncp, outvp, coord,
				invp->len);
		goto done;
	}
	/* else */

	status = NCvnrecs(outncp, nrecs);
	if(status != NC_NOERR)
		goto done;

	for( /*NADA*/; *coord < nrecs; (*coord)++)
	{
		status = NCxvarcpy(inncp, invp, coord,
				outncp, outvp, coord,
				invp->len);
		if(status != NC_NOERR)
			break;
	}
done:
	FREE_ONSTACK(coord);
	}
	return status;
}

/* no longer deprecated, used to support the 2.x interface  and also the netcdf-4 api. */
int
nc_get_att(int ncid, int varid, const char *name, void *value)
{
	int status;
	nc_type atttype;

	status = nc_inq_atttype(ncid, varid, name, &atttype);
	if(status != NC_NOERR)
		return status;

	switch (atttype) {
	case NC_BYTE:
		return nc_get_att_schar(ncid, varid, name,
			(schar *)value);
	case NC_CHAR:
		return nc_get_att_text(ncid, varid, name,
			(char *)value);
	case NC_SHORT:
		return nc_get_att_short(ncid, varid, name,
			(short *)value);
	case NC_INT:
#if (SIZEOF_INT >= X_SIZEOF_INT)
		return nc_get_att_int(ncid, varid, name,
			(int *)value);
#elif SIZEOF_LONG == X_SIZEOF_INT
		return nc_get_att_long(ncid, varid, name,
			(long *)value);
#endif
	case NC_FLOAT:
		return nc_get_att_float(ncid, varid, name,
			(float *)value);
	case NC_DOUBLE:
		return nc_get_att_double(ncid, varid, name,
			(double *)value);
	}
	return NC_EBADTYPE;
}


int
nc_put_att(
	int ncid,
	int varid,
	const char *name,
	nc_type type,
	size_t nelems,
	const void *value)
{
	switch (type) {
	case NC_BYTE:
		return nc_put_att_schar(ncid, varid, name, type, nelems,
			(schar *)value);
	case NC_CHAR:
		return nc_put_att_text(ncid, varid, name, nelems,
			(char *)value);
	case NC_SHORT:
		return nc_put_att_short(ncid, varid, name, type, nelems,
			(short *)value);
	case NC_INT:
#if (SIZEOF_INT >= X_SIZEOF_INT)
		return nc_put_att_int(ncid, varid, name, type, nelems,
			(int *)value);
#elif SIZEOF_LONG == X_SIZEOF_INT
		return nc_put_att_long(ncid, varid, name, type, nelems,
			(long *)value);
#endif
	case NC_FLOAT:
		return nc_put_att_float(ncid, varid, name, type, nelems,
			(float *)value);
	case NC_DOUBLE:
		return nc_put_att_double(ncid, varid, name, type, nelems,
			(double *)value);
	}
	return NC_EBADTYPE;
}


int
nc_get_var1(int ncid, int varid, const size_t *coord, void *value)
{
	int status;
	nc_type vartype;

	status = nc_inq_vartype(ncid, varid, &vartype); 
	if(status != NC_NOERR)
		return status;

	switch(vartype){
	case NC_CHAR:
		return nc_get_var1_text(ncid, varid, coord,
			(char *) value);
	case NC_BYTE:
		return nc_get_var1_schar(ncid, varid, coord,
			(schar *) value);
	case NC_SHORT:
		return nc_get_var1_short(ncid, varid, coord,
			(short *) value);
	case NC_INT:
		return nc_get_var1_int(ncid, varid, coord,
			(int *) value);
	case NC_FLOAT:
		return nc_get_var1_float(ncid, varid, coord,
			(float *) value);
	case NC_DOUBLE: 
		return nc_get_var1_double(ncid, varid, coord,
			(double *) value);
	}
	return NC_EBADTYPE;
}


int
nc_put_var1(int ncid, int varid, const size_t *coord, const void *value)
{
	int status;
	nc_type vartype;

	status = nc_inq_vartype(ncid, varid, &vartype); 
	if(status != NC_NOERR)
		return status;

	switch(vartype){
	case NC_CHAR:
		return nc_put_var1_text(ncid, varid, coord,
			(const char *) value);
	case NC_BYTE:
		return nc_put_var1_schar(ncid, varid, coord,
			(const schar *) value);
	case NC_SHORT:
		return nc_put_var1_short(ncid, varid, coord,
			(const short *) value);
	case NC_INT:
		return nc_put_var1_int(ncid, varid, coord,
			(const int *) value);
	case NC_FLOAT:
		return nc_put_var1_float(ncid, varid, coord,
			(const float *) value);
	case NC_DOUBLE: 
		return nc_put_var1_double(ncid, varid, coord,
			(const double *) value);
	}
	return NC_EBADTYPE;
}


int
nc_get_vara(int ncid, int varid,
	 const size_t *start, const size_t *edges, void *value)
{
	int status;
	nc_type vartype;

	status = nc_inq_vartype(ncid, varid, &vartype); 
	if(status != NC_NOERR)
		return status;

	switch(vartype){
	case NC_CHAR:
		return nc_get_vara_text(ncid, varid, start, edges,
			(char *) value);
	case NC_BYTE:
		return nc_get_vara_schar(ncid, varid, start, edges,
			(schar *) value);
	case NC_SHORT:
		return nc_get_vara_short(ncid, varid, start, edges,
			(short *) value);
	case NC_INT:
#if (SIZEOF_INT >= X_SIZEOF_INT)
		return nc_get_vara_int(ncid, varid, start, edges,
			(int *) value);
#elif SIZEOF_LONG == X_SIZEOF_INT
		return nc_get_vara_long(ncid, varid, start, edges,
			(long *) value);
#else
#error "nc_get_vara implementation"
#endif
	case NC_FLOAT:
		return nc_get_vara_float(ncid, varid, start, edges,
			(float *) value);
	case NC_DOUBLE: 
		return nc_get_vara_double(ncid, varid, start, edges,
			(double *) value);
	}
	return NC_EBADTYPE;
}

int
nc_put_vara(int ncid, int varid,
	 const size_t *start, const size_t *edges, const void *value)
{
	int status;
	nc_type vartype;

	status = nc_inq_vartype(ncid, varid, &vartype); 
	if(status != NC_NOERR)
		return status;

	switch(vartype){
	case NC_CHAR:
		return nc_put_vara_text(ncid, varid, start, edges,
			(const char *) value);
	case NC_BYTE:
		return nc_put_vara_schar(ncid, varid, start, edges,
			(const schar *) value);
	case NC_SHORT:
		return nc_put_vara_short(ncid, varid, start, edges,
			(const short *) value);
	case NC_INT:
		return nc_put_vara_int(ncid, varid, start, edges,
			(const int *) value);
	case NC_FLOAT:
		return nc_put_vara_float(ncid, varid, start, edges,
			(const float *) value);
	case NC_DOUBLE: 
		return nc_put_vara_double(ncid, varid, start, edges,
			(const double *) value);
	}
	return NC_EBADTYPE;
}

int
nc_get_varm (
	int ncid,
	int varid,
	const size_t * start,
	const size_t * edges,
	const ptrdiff_t * stride,
	const ptrdiff_t * imapp,
	void *value)
{
	int status;
	nc_type vartype;
	int varndims;
	ptrdiff_t *cvtmap = NULL;

	status = nc_inq_vartype(ncid, varid, &vartype); 
	if(status != NC_NOERR)
		return status;

	status = nc_inq_varndims(ncid, varid, &varndims); 
	if(status != NC_NOERR)
		return status;

	if(imapp != NULL && varndims != 0)
	{
		/*
		 * convert map units from bytes to units of sizeof(type)
		 */
		size_t ii;
		const ptrdiff_t szof = (ptrdiff_t) nctypelen(vartype);
		cvtmap = (ptrdiff_t *)calloc(varndims, sizeof(ptrdiff_t));
		if(cvtmap == NULL)
			return NC_ENOMEM;
		for(ii = 0; ii < varndims; ii++)
		{
			if(imapp[ii] % szof != 0)	
			{
				free(cvtmap);
				return NC_EINVAL;
			}
			cvtmap[ii] = imapp[ii] / szof;
		}
		imapp = cvtmap;
	}

	switch(vartype){
	case NC_CHAR:
		status =  nc_get_varm_text(ncid, varid, start, edges,
			stride, imapp,
			(char *) value);
		break;
	case NC_BYTE:
		status = nc_get_varm_schar(ncid, varid, start, edges,
			stride, imapp,
			(schar *) value);
		break;
	case NC_SHORT:
		status = nc_get_varm_short(ncid, varid, start, edges,
			stride, imapp,
			(short *) value);
		break;
	case NC_INT:
#if (SIZEOF_INT >= X_SIZEOF_INT)
		status = nc_get_varm_int(ncid, varid, start, edges,
			stride, imapp,
			(int *) value);
#elif SIZEOF_LONG == X_SIZEOF_INT
		status = nc_get_varm_long(ncid, varid, start, edges,
			stride, imapp,
			(long *) value);
#else
#error "nc_get_varm implementation"
#endif
		break;
	case NC_FLOAT:
		status = nc_get_varm_float(ncid, varid, start, edges,
			stride, imapp,
			(float *) value);
		break;
	case NC_DOUBLE: 
		status = nc_get_varm_double(ncid, varid, start, edges,
			stride, imapp,
			(double *) value);
		break;
	default:
		status = NC_EBADTYPE;
		break;
	}

	if(cvtmap != NULL)
	{
		free(cvtmap);
	}
	return status;
}


int
nc_put_varm (
	int ncid,
	int varid,
	const size_t * start,
	const size_t * edges,
	const ptrdiff_t * stride,
	const ptrdiff_t * imapp,
	const void *value)
{
	int status;
	nc_type vartype;
	int varndims;
	ptrdiff_t *cvtmap = NULL;

	status = nc_inq_vartype(ncid, varid, &vartype); 
	if(status != NC_NOERR)
		return status;

	status = nc_inq_varndims(ncid, varid, &varndims); 
	if(status != NC_NOERR)
		return status;

	if(imapp != NULL && varndims != 0)
	{
		/*
		 * convert map units from bytes to units of sizeof(type)
		 */
		size_t ii;
		const ptrdiff_t szof = (ptrdiff_t) nctypelen(vartype);
		cvtmap = (ptrdiff_t *)calloc(varndims, sizeof(ptrdiff_t));
		if(cvtmap == NULL)
			return NC_ENOMEM;
		for(ii = 0; ii < varndims; ii++)
		{
			if(imapp[ii] % szof != 0)	
			{
				free(cvtmap);
				return NC_EINVAL;
			}
			cvtmap[ii] = imapp[ii] / szof;
		}
		imapp = cvtmap;
	}

	switch(vartype){
	case NC_CHAR:
		status =  nc_put_varm_text(ncid, varid, start, edges,
			stride, imapp,
			(const char *) value);
		break;
	case NC_BYTE:
		status = nc_put_varm_schar(ncid, varid, start, edges,
			stride, imapp,
			(const schar *) value);
		break;
	case NC_SHORT:
		status = nc_put_varm_short(ncid, varid, start, edges,
			stride, imapp,
			(const short *) value);
		break;
	case NC_INT:
#if (SIZEOF_INT >= X_SIZEOF_INT)
		status = nc_put_varm_int(ncid, varid, start, edges,
			stride, imapp,
			(const int *) value);
#elif SIZEOF_LONG == X_SIZEOF_INT
		status = nc_put_varm_long(ncid, varid, start, edges,
			stride, imapp,
			(const long *) value);
#else
#error "nc_put_varm implementation"
#endif
		break;
	case NC_FLOAT:
		status = nc_put_varm_float(ncid, varid, start, edges,
			stride, imapp,
			(const float *) value);
		break;
	case NC_DOUBLE: 
		status = nc_put_varm_double(ncid, varid, start, edges,
			stride, imapp,
			(const double *) value);
		break;
	default:
		status = NC_EBADTYPE;
		break;
	}

	if(cvtmap != NULL)
	{
		free(cvtmap);
	}
	return status;
}

int
nc_get_vars (
	int ncid,
	int varid,
	const size_t * start,
	const size_t * edges,
	const ptrdiff_t * stride,
	void *value)
{
	return nc_get_varm (ncid, varid, start, edges,
			 stride, 0, value);
}

int
nc_put_vars (
	int ncid,
	int varid,
	const size_t * start,
	const size_t * edges,
	const ptrdiff_t * stride,
	const void *value)
{
	return nc_put_varm (ncid, varid, start, edges,
			 stride, 0, value);
}



syntax highlighted by Code2HTML, v. 0.9.1