def evaluate_pce(pc_model, pc_coeffs, germ_samples): """ Evaluate PCE at a set of samples of the germ of this PCE Input: pc_model: PC object with into about PCE pc_coeffs: numpy array with PC coefficients of the RVs to be evaluated. Each column corresponds to one RV. germ_samples: numpy array with samples of the PCE grem at which the RVs are to be evaluated. Each line is one sample. The number of colums is the number of RVs. Output: Numpy array with PCE evaluations """ # Get data set dimensions etc. n_test_samples = germ_samples.shape[0] ndim = germ_samples.shape[1] npce = pc_model.GetNumberPCTerms() # Put PC germ samples in a UQTk array std_samples_uqtk = uqtkarray.dblArray2D(n_test_samples, ndim) std_samples_uqtk.setnpdblArray(np.asfortranarray(germ_samples)) # Numpy array to store all RVs evaluated from sampled PCEs rvs_sampled = np.zeros((n_test_samples, ndim)) # Evaluate PCE for RVs in each dimension for idim in range(ndim): # Create and fill UQTk array for PC coefficients c_k_1d_uqtk = uqtkarray.dblArray1D(npce, 0.0) for ip in range(npce): c_k_1d_uqtk[ip] = pc_coeffs[ip, idim] # Create UQTk array to store outputs in rv_from_pce_uqtk = uqtkarray.dblArray1D(n_test_samples, 0.0) # Evaluate the PCEs for reach input RV at those random samples pc_model.EvalPCAtCustPoints(rv_from_pce_uqtk, std_samples_uqtk, c_k_1d_uqtk) # Put evaluated samples in full 2D numpy array for isamp in range(n_test_samples): rvs_sampled[isamp, idim] = rv_from_pce_uqtk[isamp] # return numpy array of PCE evaluations return rvs_sampled
def generate_qw(ndim, param, sp='full', type='LU'): # get quad points and weights x = uqtkarray.dblArray2D() w = uqtkarray.dblArray1D() #print 'Create an instance of Quad class' q = uqtkquad.Quad(type, sp, ndim, param) #print 'Now set and get the quadrature rule...' q.SetRule() q.GetRule(x, w) # print out x and w #print 'Displaying the quadrature points and weights:\n' #print x #print w n = len(x) # get quad points x_np = zeros((n, ndim)) x.getnpdblArray(x_np) # get quad weights w_np = zeros(n) w.getnpdblArray(w_np) xpts = array((x_np)) return xpts, w_np
def GalerkinProjection(pc_model,f_evaluations): """ Obtain PC coefficients by Galerkin Projection Input: pc_model : PC object with info about basis to project on f_evaluations: 1D numpy array (vector) with function to be projected, evaluated at the quadrature points Output: Numpy array with PC coefficients """ # Get parameters if len(f_evaluations.shape) > 1: print("This function can only project single variables for now") exit(1) # Number of quadrature points npce = pc_model.GetNumberPCTerms() nqp = f_evaluations.shape[0] # UQTk array for PC coefficients for one variable c_k_1d_uqtk = uqtkarray.dblArray1D(npce,0.0) # UQTk array for function evaluations at quadrature points for that variable f_uqtk = uqtkarray.dblArray1D(nqp,0.0) for ipt in range(nqp): f_uqtk[ipt]=f_evaluations[ipt] # Galerkin Projection pc_model.GalerkProjection(f_uqtk,c_k_1d_uqtk) # Put PC coefficients in numpy array c_k = np.zeros(npce) for ip in range(npce): c_k[ip] = c_k_1d_uqtk[ip] # Return numpy array of PC coefficients return c_k
def Unif2ArbSamples(rv_arb_in, rv_unif_in, verbose=0): """Given a set of original RV samples (rv_samples_in) from an arbitrary distribution, build a Rosenblatt map and use it to convert a set of uniformly distributed input samples (rv_samples_in) into a set of samples distributed according to the original set of RV samples. Input: rv_arb_in : numpy array with input RV samples from arbitrary distribution. Each line is a sample. The columns represent the dimensions of the input RV. rv_unif_in : numpy array with uniformly distributed samples (on [0,1]) to be mapped to the distribution of rv_arb_in. Each line is a sample. The columns represent the dimensions of the input RV and should match the ones in rv_arb_in. verbose : verbosity level (more output for higher values). Returns : numpy array with mapped samples. Same dimensionality as rv_unif_in. NOTE: * Maybe rename to InverseRosenblatt * Do a test that the input std samples are between 0 and 1 (is their std dev also 1 ??) """ # Algorithm parameters bw = -1 # KDE bandwidth for Rosenblatt (on interval 0.1) iiout = 500 # interval for output to screen # Dimensionality and number of samples of input RVs ndim_arb_in = rv_arb_in.shape[1] nsamp_arb_in = rv_arb_in.shape[0] # Set up transpose of input data for the inverse Rosenblatt transformation in a UQTk array # Rosenblatt functions want input with dimensions of [ndim,nsamp] rv_arb_in_uqtk = uqtkarray.numpy2uqtk(rv_arb_in.T) # Dimensionality and number of uniform samples to be converted to the KL RV # distribution using the inverse Rosenblatt map ndim_in = rv_unif_in.shape[1] nsamp_in = rv_unif_in.shape[0] # check that dimensionality matches the KL RV dimensionality if ndim_in != ndim_arb_in: print "Unif2ArbSamples: dimensionality of input uniform RVs does not match dimensionality of target RVs!" quit(1) # Set up numpy array for mapped sample points mappedRVs_np = np.zeros((nsamp_in, ndim_in)) # Map all samples in the uniform input to the distribution given by the data # using the inverse Rosenblatt transformation for ipt in range(nsamp_in): # print "Converting sample point #",ipt # Set up working arrays unif_uqtk_1s = uqtkarray.dblArray1D(ndim_in, 0.0) mappedRV_uqtk_1s = uqtkarray.dblArray1D(ndim_in, 0.0) # Store current sample on [0,1] in work array for idim in range(ndim_in): unif_uqtk_1s[idim] = rv_unif_in[ipt, idim] # Map each point from uniform[0,1] to the distribution given by the original samples via inverse Rosenblatt if bw > 0: uqtktools.invRos(unif_uqtk_1s, rv_arb_in_uqtk, mappedRV_uqtk_1s, bw) else: uqtktools.invRos(unif_uqtk_1s, rv_arb_in_uqtk, mappedRV_uqtk_1s) # Store results for idim in range(ndim_in): mappedRVs_np[ipt, idim] = mappedRV_uqtk_1s[idim] # Screen diagnostic output if verbose > 1: if ((ipt + 1) % iiout == 0) or ipt == 0 or (ipt + 1) == nsamp_in: print "Inverse Rosenblatt of uniform input samples:", ( ipt + 1), "/", nsamp_in, "=", (ipt + 1) * 100 / nsamp_in, "% completed" # Return samples mapped to target distribution return mappedRVs_np
def map2pce(pc_model, rvs_in, verbose=0): """Obtain PC representation for the random variables that are described by samples. Employ a Rosenblatt transformation to build a map between the input RVs and the space of the PC germ. Input: pc_model: object with properties of the PCE to be constructed rvs_in : numpy array with input RV samples. Each line is a sample. The columns represent the dimensions of the input RV. verbose : verbosity level (more output for higher values) Output: Numpy array with PC coefficients for each RV in the original rvs_in input """ # Dimensionality and number of samples of input RVs ndim = rvs_in.shape[1] nsamp = rvs_in.shape[0] # Algorithm parameters bw = -1 # KDE bandwidth for Rosenblatt (on interval 0.1) iiout = 500 # interval for output to screen # Number of PCE terms npce = pc_model.GetNumberPCTerms() # Get the default quadrature points qdpts = uqtkarray.dblArray2D() pc_model.GetQuadPoints(qdpts) totquat = pc_model.GetNQuadPoints() print "Total number of quadrature points =", totquat # Set up transpose of input data for the inverse Rosenblatt transformation in a UQTk array ydata_t = uqtkarray.dblArray2D(ndim, nsamp) ydata_t.setnpdblArray(np.asfortranarray(rvs_in.T)) # Set up numpy array for mapped quadrature points invRosData = np.zeros((totquat, ndim)) # Map all quadrature points in chosen PC set to the distribution given by the data # using the inverse Rosenblatt transformation for ipt in range(totquat): # print "Converting quadrature point #",ipt # Set up working arrays quadunif = uqtkarray.dblArray1D(ndim, 0.0) invRosData_1s = uqtkarray.dblArray1D(ndim, 0.0) # First map each point to uniform[0,1] # PCtoPC maps to [-1,1], which then gets remapped to [0,1] for idim in range(ndim): quadunif[idim] = (uqtktools.PCtoPC( qdpts[ipt, idim], pc_model.GetPCType(), pc_model.GetAlpha(), pc_model.GetBeta(), "LU", 0.0, 0.0) + 1.0) / 2.0 # Map each point from uniform[0,1] to the distribution given by the original samples via inverse Rosenblatt if bw > 0: uqtktools.invRos(quadunif, ydata_t, invRosData_1s, bw) else: uqtktools.invRos(quadunif, ydata_t, invRosData_1s) # Store results for idim in range(ndim): invRosData[ipt, idim] = invRosData_1s[idim] # Screen diagnostic output if ((ipt + 1) % iiout == 0) or ipt == 0 or (ipt + 1) == totquat: print "Inverse Rosenblatt for Galerkin projection:", ( ipt + 1), "/", totquat, "=", (ipt + 1) * 100 / totquat, "% completed" # Get PC coefficients by Galerkin projection # Set up numpy array for PC coefficients (one column for each transformed random variable) c_k = np.zeros((npce, ndim)) # Project each random variable one by one for idim in range(ndim): # UQTk array for PC coefficients for one variable c_k_1d = uqtkarray.dblArray1D(npce, 0.0) # UQTk array for map evaluations at quadrature points for that variable invRosData_1d = uqtkarray.dblArray1D(totquat, 0.0) # invRosData_1d.setnpdblArray(np.asfortranarray(invRosData[:,idim]) for ipt in range(totquat): invRosData_1d[ipt] = invRosData[ipt, idim] # Galerkin Projection pc_model.GalerkProjection(invRosData_1d, c_k_1d) # Put coefficients in full array for ip in range(npce): c_k[ip, idim] = c_k_1d[ip] # Return numpy array of PC coefficients return c_k
def Gauss2KLRVs(kl_rvs_in, gauss_rvs_in, verbose=0): """Given a set of original KL RV samples (kl_rvs_in), build a Rosenblatt map and use it to convert a set of normally distributed input samples (gauss_rvs_in) into a set of samples distributed according to the original set of KL RV samples. Input: kl_rvs_in : numpy array with input RV samples. Each line is a sample. The columns represent the dimensions of the input RV gauss_rvs_in: numpy array with normally distributed samples to be mapped to the distribution of kl_rvs_in. Each line is a sample. The columns represent the dimensions of the input RV. Number of columns should match the ones in kl_rvs_in verbose : verbosity level (more output for higher values) Returns : numpy array with mapped samples. Same dimensionality as gauss_rvs_in """ # Algorithm parameters bw = -1 # KDE bandwidth for Rosenblatt (on interval 0.1) iiout = 500 # interval for output to screen # Dimensionality and number of samples of input RVs ndim_kl = kl_rvs_in.shape[1] nsamp_kl = kl_rvs_in.shape[0] # Set up transpose of input data for the inverse Rosenblatt transformation in a UQTk array # Rosenblatt functions want input with dimensions of [ndim,nsamp] kl_samp_uqtk = uqtkarray.numpy2uqtk(kl_rvs_in.T) # Dimensionality and number of Gaussian samples to be converted to the KL RV # distribution using the inverse Rosenblatt map ndim_in = gauss_rvs_in.shape[1] nsamp_in = gauss_rvs_in.shape[0] # check that dimensionality matches the KL RV dimensionality if ndim_in != ndim_kl: print "Gauss2KLRVs: dimensionality of input Gaussian RVs does not match KL RV target dimensionality" quit(1) # Set up numpy array for mapped sample points mappedRVs_np = np.zeros((nsamp_in, ndim_in)) # Map all samples in Gaussian input to the distribution given by the data # using the inverse Rosenblatt transformation for ipt in range(nsamp_in): # print "Converting sample point #",ipt # Set up working arrays unif_uqtk_1s = uqtkarray.dblArray1D(ndim_in, 0.0) mappedRV_uqtk_1s = uqtkarray.dblArray1D(ndim_in, 0.0) # First map each point to uniform[0,1] # PCtoPC maps to [-1,1], which then gets remapped to [0,1] for idim in range(ndim_in): unif_uqtk_1s[idim] = (uqtktools.PCtoPC( gauss_rvs_in[ipt, idim], "HG", 0.0, 0.0, "LU", 0.0, 0.0) + 1.0) / 2.0 # Map each point from uniform[0,1] to the distribution given by the original samples via inverse Rosenblatt if bw > 0: uqtktools.invRos(unif_uqtk_1s, kl_samp_uqtk, mappedRV_uqtk_1s, bw) else: uqtktools.invRos(unif_uqtk_1s, kl_samp_uqtk, mappedRV_uqtk_1s) # Store results for idim in range(ndim_in): mappedRVs_np[ipt, idim] = mappedRV_uqtk_1s[idim] # Screen diagnostic output if verbose > 1: if ((ipt + 1) % iiout == 0) or ipt == 0 or (ipt + 1) == nsamp_in: print "Inverse Rosenblatt of Gaussian input samples:", ( ipt + 1), "/", nsamp_in, "=", (ipt + 1) * 100 / nsamp_in, "% completed" # Return samples mapped to target distribution return mappedRVs_np
def mkSetDbl1D(v_np): v_uqtk = uqtkarray.dblArray1D(v_np.shape[0]) v_uqtk.setnpdblArray(v_np) return v_uqtk
def mkBlnkDbl1D(): return uqtkarray.dblArray1D()
def mkSzDbl1D(dim1): return uqtkarray.dblArray1D(dim1)