#!/usr/bin/FEWXPYTHON

#######  README 
# Routine is based of Koch et al 1983 Journal of Climate and Applied MEteorology "An Interactive Barnes Objective Map Analysis Scheme for Use with Satellite and Conventional Data"
#
# When taking unevenly spaced data over a domain, the value at a grid point in the domain is estimaged as follows:
# 	Final Val = First Pass Val + Second Pass Val  (where first pass is its own array and second pass is its own separate array)
#
# Only 2 passes are needed when the proper Gamma value is used for the data set
#
# The First pass value at a grid is simply the distance weighted average of all data within 
# radius of influence - (roi) -  max distance a data pt can reside from a grid point to be factored into the estimated value at the grid point -- beyond this radius data pts are ignored
#
# The second pass value is the difference of the first pass val AT THE SOURCE DATA LOCATIONS from the actual data val -- then spread back to the grid using distance weighting.  Two ways to calc diff at src data location:
# -- a. Per Koch et al (original method) - the difference is that between PASS 1 and the src data value at every src data location 
# -- b. Per Daniels (Timothy in FE) - src data is directly distance weighted to all other src data loctions (not grid locations)- a "grid-ignorant approach".  The diff between orig val and this value is then spread over grid.
# 
# The final val is the addition of the two passes
#
# NOTE - the Koch method of the second pass differs slightly than described in the paper since this code is not using a bilinear back-interpolation (which the Koch paper does NOT imply is "THE" only method -- p 1489 last para).
#
################### Functions in order of appearance
def EXITING(t0):
	tf=datetime.datetime.now()
	dt=tf-t0
	print(str(tf) + " - exiting "+sys.argv[0]+" script - total run time "+str(dt)+" (h:m:s)")
	print("")
	print("Acceptable command line inputs: ",validsettings) 
	print("")
	exit()

def READDATASET(datafile,subdomain):
	if subdomain == "": 
		minlon=-999; maxlon=999
		minlat=-999; maxlat=999
	else:
		minlon=subdomain[0]; maxlon=subdomain[1];minlat=subdomain[2];maxlat=subdomain[3]
		if minlon >= maxlon or maxlon <= minlon or minlat >= maxlat or maxlat <= minlat:
			print ("Adjust your userdatadomain format and values: needs to be minlon,maxlon,minlat,maxlat where the mins have to be less than the maxs")
			tnow=datetime.datetime.now()
			EXITING(tnow)
		if minlon < -180 or maxlon > 360 or minlat < -90 or maxlat > 90: 
			print ("Adjust your userdatadomain values: Lats have to range from -90 to 90 and lons -180 to 360")
			tnow=datetime.datetime.now()
			EXITING(tnow)
	
	dataval=np.array([]);datalat=np.array([]);datalon=np.array([]);curline=1
	with open(datafile,"r") as FILE:
		for line in FILE:
			# remove any accidental leading/trailing white spaces
			lat=line.split(",")[0].strip()
			lon=line.split(",")[1].strip()
			val=line.split(",")[2].strip()

			# if value has an E or M take the value after that char
			if val.count("E") > 0: val=val.split("E")[1]
			if val.count("M") > 0: val=val.split("M")[1]

			# if any of these strings lat lon or val have a char thats not associated with a real number then its invalid (validdata = 0) -- example any headerline will be ignored
			validdata=bool(re.match(r"^-?\d*\.?\d*$",lat)) * bool(re.match(r"^-?\d*\.?\d*$",lon)) * bool(re.match(r"^-?\d*\.?\d*$",val))

			# finally make sure we dont have nulls (which would pass the above screening test) AND we are within domain -- if we have good data then append them in the data array as floats
			if lat != "" and lon != "" and val != "" and validdata == 1:
				if float(lon) >= minlon and float(lon) <= maxlon and float(lat) >= minlat and float(lat) <= maxlat:  
					dataval=np.append(dataval,float(val))
					datalon=np.append(datalon,float(lon))
					datalat=np.append(datalat,float(lat))

			curline +=1
	
	if len(dataval) == 0: 
		print ("There was no data read in "+datafile)
		tnow=datetime.datetime.now()
		EXITING(tnow)

	return(dataval,datalat,datalon)



def barnescore(lonsrc,latsrc,valsrc,D1londest,D1latdest,userminnumsrcptsrequired,userroi,userkappa,usergamma,usermissval):
	# need to make this a generic function because although traditionally the "source" is a set of data points and the "destination" is a set of grid points to store results during pass 1 and 2, 
	# however between pass 1 and 2 there's a step needed to back-interpolate pass1 grid back to the location of the data points -- so the grid becomes the source and the datapoint locs become the destination 
	#"src"  --> source points       (source data - data being drawn from to conduct the barnes analysis) 
	#"dest" --> destination points  (where barnes results will be stored)                   
	
	# initialize 1d array of final values -- just a placeholder to get the array shape right
	D1valdest=D1londest*0.0 

	# initialize the "while" counters	
	destindex=0; destindexmax=len(D1valdest)

	# start loop through destination indices to calculate corresponding final values to be returned
	while (destindex < destindexmax):
	     	# initialize the finalval  
		finalvaldest=usermissval

		# get the lat lon that corresponds to this index
		londest=D1londest[destindex]; latdest=D1latdest[destindex]

		# data pt screening - at the curr destination lat lon, create list of source data indices residing within x-y roi (not a true roi for now) -- beyond this x-y distance assume data too far to use
		nearbydataindices=np.where( (lonsrc >= (londest-userroi)) & (lonsrc <= (londest+userroi)) & (latsrc >= (latdest-userroi)) & (latsrc <= (latdest+userroi)) )[0]

		# if sufficient number of source data pts nearby, loop through the source data pts that survived the initial screening to see if they are within the roi       
		if len(nearbydataindices) >= userminnumsrcptsrequired:
			
			# initialize some counters
			numqualifyingpoints=0; sumofweightedvalues=0; sumofwms=0

			# for any data pt within roi, use it in the barnes analysis 
			for srcindex in nearbydataindices:
				# get lat lon src data & calc precise dist between dest & src -- above screening was initial rough cut without precision - NOTE use manual haversine calc vs module -- 30% faster!
				lonsrcpt=lonsrc[srcindex];latsrcpt=latsrc[srcindex]

				##-------------------------testing/validation which distance method accurate and fastest ---- in line calc the most accurate and fastested not the import haversine, Unit  module
				# testing haversine vs euclidian method -- there is detrimental difference using euclidian -- some some cases 10-20km different within a 1x1 deg domain east us (-76,-75,40,41)
				# distanceindegrees_euclidian=((londest-lonsrcpt)**2 + (latdest-latsrcpt)**2)**.5; print( (distanceindegrees - distanceindegrees_euclidian)/(degPerMeter*1000))
				#
				# haversine test data a=.0075356, b=.77827, c=.65793, d=.40743, e= .21616, f=55.411 degrees per https://en.wikipedia.org/wiki/Haversine_formula 2/3 down the page
				#latdest=38.898;londest=-77.037;latsrcpt=48.858;lonsrcpt=2.294;a=(1-math.cos(DTR*abs(latdest-latsrcpt)))/2;b=math.cos(DTR*latdest);c=math.cos(DTR*latsrcpt);d=(1-math.cos(DTR*abs(londest-lonsrcpt)))/2
				#e=a+(b*c*d);print(e);f=math.acos(1-2*e); print(f/DTR)
				#                  math.acos(  1 - 2 * (                     a                       + (         b             *       c                *                       d                     ) ) )/DTR  	
				#-----------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------

				distanceindegrees= math.acos(  1 - 2 * ( ((1-math.cos(DTR*abs(latdest-latsrcpt)))/2) + ( math.cos(DTR*latdest) * math.cos(DTR*latsrcpt) * ((1-math.cos(DTR*abs(londest-lonsrcpt)))/2) ) ) )/DTR
				r2=distanceindegrees**2 # distance squared needed in weighting per paper's equations page 1489 (eq 1 and eq 7)

				# if the precise distance is within roi, let the data pt contribute to the barnes analysis and increment the tally of qualifying points
				if distanceindegrees <= userroi :
					valatpt=valsrc[srcindex]
					numqualifyingpoints += 1
					wm=math.exp(-r2/(userkappa*usergamma)) # per Koch paper - the weight using same "wm" and "r2" nomenclature (function of distance and kappa)
					weightedvalue=valatpt*wm; sumofweightedvalues=sumofweightedvalues+weightedvalue; sumofwms=sumofwms+wm

			# Now that the loop is completed, if we have min num data pts required that fell within search rad, update D1valdest at lat lon, otherwise it remains stored with nan
			if numqualifyingpoints >= userminnumsrcptsrequired: finalvaldest=sumofweightedvalues/sumofwms
		
			D1valdest[destindex]=finalvaldest

		destindex += 1

	return(D1valdest)



############## import modules
import numpy as np					# for array math
import math						# for use of pi, cosine, and populating vars with nan
import re						# for scanning strings to ensure we dont treat strings with letters as numbers
import xarray as xr					# for consolidating val, lats, lons into a data set and writing to grib2
import os               				# for access to system env vars
import sys						# to handle command line input
import datetime						# for time elapsed tracking 
import subprocess                                       # so we can run system commands
from cfgrib.xarray_to_grib import to_grib               # to write out results to a grib file



############## set constants
t0=datetime.datetime.now(); print (""); print(str(t0)+" - "+sys.argv[0]+" script starting (h:m:s)"); tlast=t0
Rearth=6378137			# radius of earth in meters
pi=math.pi			# pi 3.1415.....
degPerMeter=360/(2*pi*Rearth)	# degrees lon per meter
DTR=math.pi/180			# degrees to radian 
missingval=np.nan		# missing value for pre-write cals - THIS IS NOT WHAT THE USER CAN REASSIGN MISSING VALUES TO -- this is the value used through array calcs -- best to leave it as np.nan
replacemiss=""			# if replacemiss gets set by user then the last step before writing out to grib2 will be to replace any missing values with this value
############# end setting constants



############ user settings 
datafile       = "" # leave this empty -- force user to define on command line
outfile        = "" # leave this empty -- force user to define on command line
userdatadomain = "" # data domain user wishes to focus analysis    #-76.0,-75.0,40.0,41.0 #-126.0,-65.0,24.0,50.0 # -85.0,-84.0,37.0,38.0 
usergriddomain = "" # domain we will conduct analysis -- should be no larger than the userdatadomain -- and small enough so some data points lie beyond grid domain
userDN         = "" # data spacing in deg 
userDX         = "" # horizontal grid resolution in deg  -- should be no LARGER than DN/2
userkappa      = "" # barnes weighting parameter (really no need to set this -- is a fn of DN) we allow it to be manually set just to allow user to see how it influences the result
userroi        = "" # distance in degrees between grid point and data point, greater than this distance real data will not influence the estimated value at the grid point
usergamma      = "" # barnes convergence parameter, real num between 0 and 1 that influences shape of the gamma distribution (smaller means tighter, less smooth and vice versa) and is applied to kappa during second pass only
usermindatapts = "" # has to be at least 1 and ideally per Koch 1983 should be no less than 3 
userbackinterp = "" # defaults to second pass backinterpolation  method used in Koch paper -- the other method is "daniels" -- see comments at top of code for details
remapres       = "" # remap output to user specified resolution -- leave empty -- default is to NOT remap
debug          = "" # 1 =yes 0= no
# read from commandline
argmax=len(sys.argv); currarg=1; validsettings=["help", "datafile", "outfile", "datadomain", "griddomain", "dn", "dx", "kappa", "roi", "gamma", "mindatapts", "backinterp", "remap", "debug", "replacemiss"]
notrecognized=""
notdefinedright=""
while currarg < argmax:
	if sys.argv[currarg].count("help"): 
		print(" ")
		print("Python Barnes Objective Analysis Routine Help")
		print("---------------------------------------------")
		print("Recommended approach- initially only provide 3 inputs")
		print("  1. Specify the input file -- must be csv in the format of lat,lon,value with a header (example datafile=/some/path/data.csv)")
		print("  2. Specify the output file (example outfile=result.grb2).  Although the name supplied will be used, the content will be grib2.")
		print("  3. Specify the data domain")
		print("  Let the routine calculate the rest of the settings including highest recommended grid res given data spacing over domain")
		print("  and review the results -- then if needed make viable adjustments -- some common options are:")
		print("     a. If domain is really small (2 deg or less), set the grid domain to data domain")
		print("     b. Adjust gamma to between .01 to .20 to increase convergence of grid point estimations to the nearest data point")
		print(" ")
		print("Tips...")
		print("  1. Setting grid domain to data domain for small data domain areas will allow analysis to more closely match what human visualizes")
		print("  2. If you desire a higher horiz res output grib2 file than what is objectively determined, the fastest way to do this is let the")
		print("     routine produce the output grib file at the objectively calculated horiz resolution -- invoke the REMAP option which will use")
		print("     wgrib2 -new_grid option to remap the result to the horiz grid res of your choice -- much much faster than forcing a grid res higher than data supports")
		print("  3. The minimum number of data points when less than 3 will lead to holes if roi <= 2")
		print("  4. Make gamma smaller than .2 to converge to a data point (.01 to .1 is good)")
		print("  5. For large domains with rich data, the Daniels method is faster but there is no clear winner for small domains")
		print(" ")
		print("Settings Description (* = typically no need to set, let routine calculate its settings)")
		print("datafile    :     standard filename INCLUDING explicit (no env vars) path.                                 Example     datafile=/some/path/filename")
		print("outfile     :     ditto -- name and extension arbitrary -- but content always grib2")
		print("datadomain  :     minlon,maxlon,minlat,maxlat where commas are literal.                                    Example   datadomain=-76.0,-75,40,41.5")
		print("griddomain* :     ditto")
		print("dn*         :     average minimum distance in degrees between data points in datadomain.                   Example           dn=1.25")
		print("dx*         :     horiz grid res in deg on CED proj (converted from meters if value ends in m)             Example           dx=.05 or dx=5000m")
		print("kappa*      :     weighting function in degrees squared.                                                   Example        kappa=.23")
		print("roi*        :     radius of influence in degrees (beyond this data will not influence grid pt value.       Example          roi=1.4")
		print("gamma*      :     convergence value to force grid pt to data pt value  (used in 2nd pass) between 0 & 1.   Example        gamma=0.3")
		print("mindatapts* :     minimum data pts required within roi for a grid pt to be assigned a value (must be > 0). Example   mindatapts=3")
		print("backinterp* :     how to generate 2nd pass diff between orig data & interp vals at data locs.              Example   backinterp=daniels (only option, default is original method)")
		print("remap*      :     horiz res to remap output to lat lon grid (converted from meters if value ends in m)     Example        remap=.25 or remap=1000m")
		print("debug*      :     flag to turn on debug info -- either 1 or 0.                                             Example        debug=1")
		print("replacemiss :     replace missing vals to entered value                                                    Example  replacemiss=0")
		print(" ")
		print("Approach based on Koch et al 1983 Sep Journal of Climate and Applied Meteorology p 1487-1503")
		print("")
		tf=datetime.datetime.now()
		EXITING(tf)

	if sys.argv[currarg].count("datafile="):         datafile=sys.argv[currarg].split("datafile=")[1]
	if sys.argv[currarg].count("outfile="):           outfile=sys.argv[currarg].split("outfile=")[1]
	if sys.argv[currarg].count("datadomain="): userdatadomain=tuple( [ float(num_str) for num_str in sys.argv[currarg].split("datadomain=")[1].split(",") ] )
	if sys.argv[currarg].count("dn="):
		dum=sys.argv[currarg].split("dn=")[1]
		if dum != "":
			userDN=float(dum)
		else:
			notdefinedright=notdefinedright+" dn"

	if sys.argv[currarg].count("dx="):	  	   userDX=sys.argv[currarg].split("dx=")[1] 	# do not convert to float yet -- do this later as we will need to check if it was entered in deg or meters
	if sys.argv[currarg].count("griddomain="): usergriddomain=tuple( [ float(num_str) for num_str in sys.argv[currarg].split("griddomain=")[1].split(",") ] )
	if sys.argv[currarg].count("kappa="):
		dum=sys.argv[currarg].split("kappa=")[1]
		if dum != "":
			userkappa=float(dum)
		else:
			notdefinedright=notdefinedright+" kappa="

	if sys.argv[currarg].count("roi="):
		dum=sys.argv[currarg].split("roi=")[1]
		if dum != "":
			userroi=float(dum)
		else:
			notdefinedright=notdefinedright+" roi="

	if sys.argv[currarg].count("gamma="):
		dum=sys.argv[currarg].split("gamma=")[1]
		if dum != "":
			usergamma=float(dum)
		else:
			notdefinedright=notdefinedright+" gamma="

	if sys.argv[currarg].count("mindatapts="): usermindatapts=int(sys.argv[currarg].split("mindatapts=")[1])
	if sys.argv[currarg].count("backinterp="): userbackinterp=sys.argv[currarg].split("backinterp=")[1].lower()
	if sys.argv[currarg].count("remap="):            remapres=sys.argv[currarg].split("remap=")[1]	# do not convert to float yet -- do this later as we will need to check if it was entered in deg or meters
	if sys.argv[currarg].count("debug="):               debug=int(sys.argv[currarg].split("debug=")[1])
	if sys.argv[currarg].count("replacemiss="): 
		dum=sys.argv[currarg].split("replacemiss=")[1]
		if dum != "": 
			replacemiss=float(dum)
		else:
			notdefinedright=notdefinedright+" replacemiss"
			

	if not sys.argv[currarg].split("=")[0] in validsettings: notrecognized=notrecognized+" "+sys.argv[currarg].split("=")[0]
	currarg += 1

# warn user and exit if an invalid or mistyped setting entered
if notrecognized != "": 
	print ("")
	print ("The following settings you input are not recognized (could be misspelled):")
	print (notrecognized)
	print ("(type "+sys.argv[0]+" help  for a list of valid settings)")
	print ("")
	tf=datetime.datetime.now()
	EXITING(tf)

if notdefinedright != "":
	print ("")
	print ("The following settings you input are set improperly:")
	print (notdefinedright)
	print ("")
	tf=datetime.datetime.now()
	EXITING(tf)
########## end user settings



############ basic error checking for file name input by user
if outfile == "":
	print ("need an output file name -- specify this on the command line using outfile=")
	print ("")
	tf=datetime.datetime.now()
	EXITING(tf)

if datafile == "":
	print ("need an input file name -- specify this on the command line using datafile=")
	print ("")
	tf=datetime.datetime.now()
	EXITING(tf)

if not os.path.exists(datafile):
	print("Input data file "+ datafile + " was not found")
	print ("")
	tf=datetime.datetime.now()
	EXITING(tf)
############end basic error checking of filenames



############ OBJECTIVE SETTINGS based on data file contents (can be overridden by user input)
tnow=datetime.datetime.now(); print(str(tnow)+ " - reading "+datafile); tlast=tnow
dataval,datalat,datalon=READDATASET(datafile,userdatadomain)
tnow=datetime.datetime.now(); dt=tnow-tlast;print(str(tnow)+ " - finished reading "+datafile+" in "+str(dt)+"(h:m:s)");tlast=tnow

# data domain based on raw data and also the area used to calc DN based on a uniform data spacing (p 1492 eq 12 Koch 1983) -- latter will be overridden later if user set their own data domain
minlon_rawdata=np.min(datalon); maxlon_rawdata=np.max(datalon)
minlat_rawdata=np.min(datalat); maxlat_rawdata=np.max(datalat)
datadomain_objective=[minlon_rawdata,maxlon_rawdata,minlat_rawdata,maxlat_rawdata]
area_rawdata=(maxlon_rawdata-minlon_rawdata)*(maxlat_rawdata-minlat_rawdata)

# avg min distance between data points which drives DN - fast calc (vs d=sqrt(dlon^2 + dlat^2) looped over each point) is tocreate arrays of the components of dlon and dlat and then np math it from there
dlons=datalon[:,np.newaxis]-datalon[np.newaxis,:]  # create array of dlon from each other - get 1 row per pt, each holding list of dlons relative to that pt - always 0 when it subtracts from itself (diagnoal elements)
dlats=datalat[:,np.newaxis]-datalat[np.newaxis,:]  # ditto with lats
distances=np.sqrt(dlons**2 + dlats**2)             # use the two arrays dlon and dlat to quickly create a distance matrix
np.fill_diagonal(distances,np.inf)                 # Since we eventually want to hone in on the min val per row and all vals are +, replace all zeros with positive infinity so zeros to give a false minimum
min_distance=np.min(distances,axis=1)              # find the minimum in each row and put that in its own array (you get a 1d array with as many elements as you had rows)
numdatapts=len(dataval)				   # the number of data points we have
DN_objective_calculated=np.sum(min_distance)/numdatapts # take the avg of these minimums (sum of the minimums divided by number of data pts) to get your avg minimum data spacing

# for ref only: DN calculated via the "random method" Koch 1983 p1492 eq 12 (dn if all data points were uniformly spread across the grid)
DN_objective_random=(area_rawdata**.5)*((1 + numdatapts**.5)/(numdatapts-1))

# calc dx (and assume CED so dy is the same) -- then shrink grid domain on each side by DN_objective_calculated/2 deg (1 dx .. which is 1 grid pt per Koch 1983 p1493 per fig 4 and implied in section d 3rd sentence)
DX_objective=DN_objective_calculated/2             # dx should be no LARGER than this -- can be smaller though --- see Koch 1983 p 1493 section d. 3rd sentence
griddomain_objective=[minlon_rawdata+DN_objective_calculated,maxlon_rawdata-DN_objective_calculated,minlat_rawdata+DN_objective_calculated,maxlat_rawdata-DN_objective_calculated]
griddomain_objective_rounded=[round(minlon_rawdata+DN_objective_calculated,5),round(maxlon_rawdata-DN_objective_calculated,5),round(minlat_rawdata+DN_objective_calculated,5),round(maxlat_rawdata-DN_objective_calculated,5)]

# rest of the settings 
minnumdataptsrequired_objective=3                            # see Koch 1983 p 1493 para 2
kappa_objective=5.0514573*(2*DN_objective_calculated/pi)**2  # see Koch 1983 p 1492 equation 13
roi_objective=(20*kappa_objective)**.5                       # see koch 1983 p 1493 second to last sentence of section e
gamma_objective=.2                                           # see koch 1983 p 1492 section c
############## end objective settings ##############################################################



#################### finalize settings - initially adopt all the objective settings then override with any user specified settings
datadomain            = datadomain_objective
DN                    =         DN_objective_calculated
DX                    =         DX_objective
griddomain            = griddomain_objective
kappa                 =      kappa_objective
roi                   =        roi_objective
gamma                 =      gamma_objective
minnumdataptsrequired = minnumdataptsrequired_objective

### Override the above is user input suppliec
if userdatadomain != "": datadomain            = userdatadomain 
if userDN         != "": DN                    = userDN
if userDX         != "": 
	DX=userDX.lower()  # put this to lower case in case there is a capital M for meters by accident from the user
	if DX.count("m"): DX=float(DX.split("m")[0])*degPerMeter # convert it to deg if they entered meters
if usergriddomain != "": griddomain            = usergriddomain
if userkappa      != "": kappa                 = userkappa
if userroi        != "": roi                   = userroi
if usergamma      != "": gamma                 = usergamma
if usermindatapts != "": minnumdataptsrequired = usermindatapts

#---if userDN is defined, recalc user settings dependant on DN that were left undefined 
if userDN != "":
	if userDX         == "": DX=DN/2
	if usergriddomain == "": griddomain=[minlon_rawdata+DN,maxlon_rawdata-DN,minlat_rawdata+DN,maxlat_rawdata-DN] 
	if userkappa      == "": kappa=5.0514573*(2*DN_objective_calculated/pi)**2
	if userroi        == "": roi=(20*kappa)**.5

# update min/maxlon and min/maxlaat and the printable version of griddomain in case 
minlon=griddomain[0];maxlon=griddomain[1];minlat=griddomain[2];maxlat=griddomain[3]
griddomain_rounded=[round(minlon,5),round(maxlon,5),round(minlat,5),round(maxlat,5)]
######################## end finalize settings



################ output finalized settings to user (include the objective calcs as reference)
printDN="same" 
printDX="same"
printdomain="same"
printkappa="same"
printroi="same"
printgamma="same"
printreqpts="same"
if                    DN != DN_objective_calculated         : printDN     = str(round(DN_objective_calculated,5))
if                    DX != DX_objective                    : printDX     = str(round(DX_objective,5))
if            griddomain != griddomain_objective            : printdomain = str(griddomain_objective_rounded)
if                 kappa != kappa_objective                 : printkappa  = str(round(kappa_objective,5))
if                   roi != roi_objective                   : printroi    = str(round(roi_objective,5))
if                 gamma != gamma_objective                 : printgamma  = str(round(roi_objective,5))
if minnumdataptsrequired != minnumdataptsrequired_objective : printreqpts = str(minnumdataptsrequired_objective)
print ("")
print ("Used settings with comparison to idealized settings dictated by contents of "+datafile+" having "+str(numdatapts)+" data points in domain "+str(datadomain)) 
print ("dn (avg data spacing in deg)            : "+str(round(DN,5))                 + " (idealized objective setting="+printDN+", idealized random method="+str(round(DN_objective_random,5))+")")
print ("dx (horiz grid spacing in deg)          : "+str(round(DX,5))                 + " (idealized objective setting="+printDX+")")
print ("grid domain                             : "+str(griddomain_rounded)          + " (idealized objective setting="+printdomain+")")
print ("kappa (in deg^2)                        : "+str(round(kappa,5))              + " (idealized objective setting="+printkappa+")")
print ("roi (in deg)                            : "+str(round(roi,5))                + " (idealized objective setting="+printroi+")")
print ("gamma pass 2 value (always 1 on pass 1) : "+str(gamma)                       + " (idealized objective setting="+printgamma+")")
print ("minimum number of datapts required      : "+str(minnumdataptsrequired)       + " (idealized objective setting="+printreqpts+")")
print ("")
#################################################### end print out summary



######## Final ERROR CHECKING of barnes settings - can only do this after settings finalized
if DN <= 0:
	print ("Check your dn setting - cant be 0 or less")
	print ("")
	tnow=datetime.datetime.now()
	EXITING(tnow)

if DX > DX_objective: print ("CAUTION -  DX shouldnt be any larger than DN/2 - you may want to adjust that setting")

if roi <= 0:
	print ("Check your roi setting - cant be 0 or less")
	print ("")
	tnow=datetime.datetime.now()
	EXITING(tnow)

if roi < DN: print ("Caution - your roi setting is less than the average spacing of your data -- means some grid points wont have data qualifying as 'sufficiently nearby' to perform the analysis")

if gamma <= 0 or gamma >= 1:
	print ("Check your gamma setting - has to be > 0 and < 1")
	print ("")
	tnow=datetime.datetime.now()
	EXITING(tnow)

if minnumdataptsrequired < 1:
	print ("Check your minnumdatapts setting - cant be less than 1")
	print ("")
	tnow=datetime.datetime.now()
	EXITING(tnow)

if userbackinterp != "" and userbackinterp != "daniels":
	print ("Check your backinterp setting - can only be set with the string 'daniels'")
	print ("")
	tnow=datetime.datetime.now()
	EXITING(tnow)
############## end final error checking of barnes settings



####### GRID creation ######################
# Create lists of lons and lats based on domain -- being sure lat list is stored from highest to lowest to visually correspond to a map if data printed on a screen during debugging
gridlonlist=np.arange(minlon,maxlon+DX,DX);         numxpts=len(gridlonlist)
gridlatlist=np.arange(maxlat,minlat-DX,-DX);        numypts=len(gridlatlist); numpts=numxpts*numypts

# create 1d arrays of gridlon,gridlat,gridvals so later looping is fast by exploiting that 3 arrays linked by same index -- first mesh lat & lon arrays into 2d to ensure num elements in grid properly set, then flatten each to 1d
lon_mesh,lat_mesh=np.meshgrid(gridlonlist,gridlatlist)
D1gridlons=lon_mesh.flatten(); D1gridlats=lat_mesh.flatten(); D1gridvals=D1gridlats*0.0

# get a print-ready version of the grid domain
griddomain_rounded=[round(minlon,5),round(maxlon,5),round(minlat,5),round(maxlat,5)]

# display resulting grid info to user
tnow=datetime.datetime.now()
print(str(tnow)+" - Grid created: "+str(round(DX,5))+" dx,dy deg ("+str(round(DX/degPerMeter,1))+"m) horiz res over domain "+str(griddomain_rounded)+" with "+str(numpts)+" pts in a "+str(numypts)+"X"+str(numxpts)+" rowXcol grid)")
##### end grid creation



#######################temporary source data debugging
if debug == 1:
	print(" ")
	print("---data vals lats lons---")
	print(dataval)
	print(datalat)
	print(datalon)
	print(" ")
	print("grid lat lons")
	print(np.round(gridlatlist,decimals=2)); print(np.round(gridlonlist,decimals=2))
	print (" ")
################ end debugging




#### Barnes passes -- because the core method to calc a distance weight is its own function, no while loops needed
tnow=datetime.datetime.now(); print(str(tnow)+ " - begining Barnes interpolation passes (h:m:s)"); tlast=tnow

# pass 1 - simple distance weighting of nearby data to each grid point (gamma = 1)
pass1=barnescore(datalon,datalat,dataval,D1gridlons,D1gridlats,minnumdataptsrequired,roi,kappa,1,missingval)
tnow=datetime.datetime.now(); dt=tnow-tlast;print(str(tnow)+ " - pass 1 completed in "+str(dt)+" (h:m:s)"); tlast=tnow

# intermediate step to back interpolate pass 1 to the locations of the source data (gamma still = 1)
if userbackinterp == "daniels": 
	methodtype="Daniels method"
	backinterp=barnescore(datalon,datalat,dataval,datalon,datalat,minnumdataptsrequired,roi,kappa,1,0)
else:
	methodtype="Koch (original) method"
	backinterp=barnescore(D1gridlons,D1gridlats,pass1,datalon,datalat,minnumdataptsrequired,roi,kappa,1,0)
tnow=datetime.datetime.now(); dt=tnow-tlast;print(str(tnow)+ " - backinterp calculation using "+methodtype+" completed in "+str(dt)+" (h:m:s)"); tlast=tnow

# pass 2 -  calculate the difference between the source data at its location from what pass1 estimated at the same loction -- and distance weight those difference across the grid (gamma = prefereed gamma val < 1)
pass2=barnescore(datalon,datalat,dataval-backinterp,D1gridlons,D1gridlats,minnumdataptsrequired,roi,kappa,gamma,missingval)
tnow=datetime.datetime.now(); dt=tnow-tlast;print(str(tnow)+ " - pass 2 completed in "+str(dt)+" (h:m:s)");tlast=tnow

# final result = pass 1 (the first guess) fine-tuned with the difference grid
finalresults=pass1+pass2

# debugging
if debug == 1:
	vispass1=np.round(pass1.reshape((numypts,numxpts)),decimals=2); print("pass1");print(vispass1)
	print("backinterp");print(np.round(backinterp,decimals=2))
	vispass2=np.round(pass2.reshape((numypts,numxpts)),decimals=2); print("pass2");print(vispass2)
	print("finalresults = pass1 + pass2");print(np.round(finalresults.reshape((numypts,numxpts)),decimals=2))
### end Barnes passes ##################



############## write routine
# first - replace nans with a user defined value if they set it -- do it here before we set up xarray for writing to grib2 file
if replacemiss != "": 
	tnow=datetime.datetime.now(); print(str(tnow)+" - nans being replaced with "+str(replacemiss)); tlast=tnow
	finalresults[np.isnan(finalresults)]=replacemiss


# put finalresults in 2d form for writing by cfgrib 
finalresults_2d=finalresults.reshape(numypts,numxpts)
if debug == 1:
	print(finalresults_2d.shape)
	print("............")
	print("---final---")
	print(finalresults_2d)
	print("---latlist---")
	print(gridlatlist)
	print(gridlatlist.shape)
	print("---lonlist---")
	print(gridlonlist)
	print(gridlonlist.shape)
xbarnes=xr.DataArray(data=finalresults_2d,dims=['latitude','longitude'],coords={'latitude':gridlatlist,'longitude':gridlonlist},) # dont need to assign anything else 
if debug == 1:
	print("---xbarnes---")
	print(xbarnes)
	print("==============")


xwrite=xbarnes.to_dataset(name="result") 				# "name" can be anything
if debug == 1:
	print("")
	print(xwrite)
	print("")
to_grib(xwrite,outfile,grib_keys={'edition':2}, no_warn=True)		# this is all that is needed -- will default to writing it as SFC TMP ANL with a date of 20070323 12-- which usr can change using wgrib2 

tnow=datetime.datetime.now()
msg2=""
if not os.path.exists(outfile): 
	msg=" did not get written"
else:
	if os.path.getsize(outfile) == 0: 
		msg = outfile + " created but was zero bytes"
	else:
		#### remap if user specified as such
		remapmsg=" "
		if remapres != "":  # user specified some kind of remap -- either in deg or km
			remapres=remapres.lower()  # put this to lower case in case there is a capital M for meters by accident from the user
			if remapres.count("m"): remapres=float(remapres.split("m")[0])*degPerMeter # convert it to deg if they entered meters
			# warn user if remap is coarser or lower in res than the original --- allow it but caution against it
			if remapres >= DX:
				print("                             WARNING - your remapped resolution "+str(round(remapres,5))+" is the same or coarser than the original "+str(round(DX,5))+" and not recommended -- proceeding anyway")
			# using wgrib2 -grid get teh bounds of the grib data in outfile
			wgriblons=re.sub(r"[\t]*","",subprocess.run(["wgrib2 -grid "+outfile+" | grep lon | grep to | awk '{print $2,$4}'"],shell=True,capture_output=True,text=True).stdout.rstrip()).split("\n")[0].split()
			wgriblats=re.sub(r"[\t]*","",subprocess.run(["wgrib2 -grid "+outfile+" | grep lat | grep to | awk '{print $2,$4}'"],shell=True,capture_output=True,text=True).stdout.rstrip()).split("\n")[0].split()
			# calc the number of new lat lons using the bounds and the new remap resolution
			numwgriblons=int((float(wgriblons[1])-float(wgriblons[0]))/remapres)+1
			numwgriblats=int((float(wgriblats[0])-float(wgriblats[1]))/remapres)+1
			# create a new target file name that indicates the original and remapped resolution in deg
			newoutfile=outfile.split(".")[0]+"_origresdeg_"+str(round(DX,5))+"_remappedtodeg_"+str(round(remapres,5))+".grb2"
			remapmsg=" (and remapped to "+str(round(remapres,5))+" horiz res via wgrib2 in "+newoutfile+") "
			# create the wgrib2 remap command
			remapcmd="wgrib2 "+outfile+" -new_grid latlon "+str(wgriblons[0])+":"+str(numwgriblons)+":"+str(remapres)+" "+str(wgriblats[1])+":"+str(numwgriblats)+":"+str(remapres)+" "+str(newoutfile)
			print(remapcmd)
			# execute the command 
			remapresult=re.sub(r"[\t]*","",subprocess.run([remapcmd],shell=True,capture_output=True,text=True).stdout.rstrip()).split("\n")
		msg  = outfile + " written by cfgrib"+remapmsg+"with its default placeholder attributes SFC TMP ANL at 2007032312"
		msg2 = "(if plotting, best to use domain no larger than "+str(datadomain_objective)+ " however data only was gridded in "+str(griddomain_rounded)+")"
print (str(tnow)+" - "+msg)
if msg2 != "": print("                             "+msg2)
############ end write routine ###########################



############## Exit routine
EXITING(t0)
