/** * Intialize lrthresh data * * @param dims_decom - dimensions with levels at LEVEL_DIMS * @param randshift - randshift boolean * @param mflags - selects which dimensions gets reshaped as the first dimension in matrix * @param blkdims - contains block dimensions for all levels * */ static struct lrthresh_data_s* lrthresh_create_data(const long dims_decom[DIMS], bool randshift, unsigned long mflags, const long blkdims[MAX_LEV][DIMS], float lambda, bool noise, int remove_mean) { PTR_ALLOC(struct lrthresh_data_s, data); SET_TYPEID(lrthresh_data_s, data); data->randshift = randshift; data->mflags = mflags; data->lambda = lambda; data->noise = noise; data->remove_mean = remove_mean; // level dimensions md_copy_dims(DIMS, data->dims_decom, dims_decom); md_calc_strides(DIMS, data->strs_lev, dims_decom, CFL_SIZE); // image dimensions data->levels = dims_decom[LEVEL_DIM]; md_select_dims(DIMS, ~LEVEL_FLAG, data->dims, dims_decom); md_calc_strides(DIMS, data->strs, data->dims, CFL_SIZE); // blkdims for(long l = 0; l < data->levels; l++) { for (long i = 0; i < DIMS; i++) data->blkdims[l][i] = blkdims[l][i]; } return PTR_PASS(data); }
static struct maps_data* maps_create_data(const long max_dims[DIMS], unsigned int sens_flags, const complex float* sens, bool gpu) { struct maps_data* data = xmalloc(sizeof(struct maps_data)); // maximal dimensions md_copy_dims(DIMS, data->max_dims, max_dims); // sensitivity dimensions md_select_dims(DIMS, sens_flags, data->mps_dims, max_dims); md_calc_strides(DIMS, data->strs_mps, data->mps_dims, CFL_SIZE); md_select_dims(DIMS, ~MAPS_FLAG, data->ksp_dims, max_dims); md_calc_strides(DIMS, data->strs_ksp, data->ksp_dims, CFL_SIZE); md_select_dims(DIMS, ~COIL_FLAG, data->img_dims, max_dims); md_calc_strides(DIMS, data->strs_img, data->img_dims, CFL_SIZE); #ifdef USE_CUDA complex float* nsens = (gpu ? md_alloc_gpu : md_alloc)(DIMS, data->mps_dims, CFL_SIZE); #else assert(!gpu); complex float* nsens = md_alloc(DIMS, data->mps_dims, CFL_SIZE); #endif md_copy(DIMS, data->mps_dims, nsens, sens, CFL_SIZE); data->sens = nsens; data->norm = NULL; return data; }
const struct linop_s* linop_zfinitediff_create(unsigned int D, const long dims[D], long diffdim, bool circular) { PTR_ALLOC(struct zfinitediff_data, data); SET_TYPEID(zfinitediff_data, data); data->D = D; data->dim_diff = diffdim; data->do_circdiff = circular; data->dims_in = *TYPE_ALLOC(long[D]); data->dims_adj = *TYPE_ALLOC(long[D]); data->strides_in = *TYPE_ALLOC(long[D]); data->strides_adj = *TYPE_ALLOC(long[D]); md_copy_dims(D, data->dims_in, dims); md_copy_dims(D, data->dims_adj, dims); md_calc_strides(D, data->strides_in, data->dims_in, CFL_SIZE); if (!data->do_circdiff) data->dims_adj[data->dim_diff] -= 1; md_calc_strides(D, data->strides_adj, data->dims_adj, CFL_SIZE); const long* dims_adj = data->dims_adj; const long* dims_in = data->dims_in; return linop_create(D, dims_adj, D, dims_in, CAST_UP(PTR_PASS(data)), zfinitediff_apply, zfinitediff_adjoint, zfinitediff_normal, NULL, zfinitediff_del); }
static struct ufft_data* ufft_create_data(const long ksp_dims[DIMS], const long pat_dims[DIMS], const complex float* pat, unsigned int flags, bool use_gpu) { PTR_ALLOC(struct ufft_data, data); SET_TYPEID(ufft_data, data); data->flags = flags; data->use_gpu = use_gpu; md_copy_dims(DIMS, data->pat_dims, pat_dims); md_copy_dims(DIMS, data->ksp_dims, ksp_dims); md_calc_strides(DIMS, data->pat_strs, pat_dims, CFL_SIZE); md_calc_strides(DIMS, data->ksp_strs, ksp_dims, CFL_SIZE); #ifdef USE_CUDA data->pat = (use_gpu ? md_alloc_gpu : md_alloc)(DIMS, data->pat_dims, CFL_SIZE); #else data->pat = md_alloc(DIMS, data->pat_dims, CFL_SIZE); #endif md_copy(DIMS, data->pat_dims, data->pat, pat, CFL_SIZE); data->fft_op = linop_fftc_create(DIMS, ksp_dims, flags); return PTR_PASS(data); }
static struct linop_s* linop_gdiag_create(unsigned int N, const long dims[N], unsigned int flags, const complex float* diag, bool rdiag) { PTR_ALLOC(struct cdiag_s, data); SET_TYPEID(cdiag_s, data); data->rmul = rdiag; data->N = N; PTR_ALLOC(long[N], dims2); PTR_ALLOC(long[N], dstrs); PTR_ALLOC(long[N], strs); long ddims[N]; md_select_dims(N, flags, ddims, dims); md_copy_dims(N, *dims2, dims); md_calc_strides(N, *strs, dims, CFL_SIZE); md_calc_strides(N, *dstrs, ddims, CFL_SIZE); data->dims = *PTR_PASS(dims2); data->strs = *PTR_PASS(strs); data->dstrs = *PTR_PASS(dstrs); data->diag = diag; // make a copy? #ifdef USE_CUDA data->gpu_diag = NULL; #endif return linop_create(N, dims, N, dims, CAST_UP(PTR_PASS(data)), cdiag_apply, cdiag_adjoint, cdiag_normal, NULL, cdiag_free); }
const struct operator_s* nufft_precond_create(const struct linop_s* nufft_op) { const auto data = CAST_DOWN(nufft_data, linop_get_data(nufft_op)); PTR_ALLOC(struct nufft_precond_data, pdata); SET_TYPEID(nufft_precond_data, pdata); assert(data->conf.toeplitz); int N = data->N; int ND = N + 1; pdata->N = N; pdata->cim_dims = *TYPE_ALLOC(long[ND]); pdata->pre_dims = *TYPE_ALLOC(long[ND]); pdata->cim_strs = *TYPE_ALLOC(long[ND]); pdata->pre_strs = *TYPE_ALLOC(long[ND]); md_copy_dims(ND, pdata->cim_dims, data->cim_dims); md_select_dims(ND, data->flags, pdata->pre_dims, pdata->cim_dims); md_calc_strides(ND, pdata->cim_strs, pdata->cim_dims, CFL_SIZE); md_calc_strides(ND, pdata->pre_strs, pdata->pre_dims, CFL_SIZE); pdata->pre = compute_precond(pdata->N, pdata->pre_dims, pdata->pre_strs, data->psf_dims, data->psf_strs, data->psf, data->linphase); pdata->fft_op = linop_fft_create(pdata->N, pdata->cim_dims, data->flags); const long* cim_dims = pdata->cim_dims; // need to dereference pdata before PTR_PASS return operator_create(N, cim_dims, N, cim_dims, CAST_UP(PTR_PASS(pdata)), nufft_precond_apply, nufft_precond_del); }
int main_repmat(int argc, char* argv[]) { mini_cmdline(argc, argv, 4, usage_str, help_str); long in_dims[DIMS]; long out_dims[DIMS]; complex float* in_data = load_cfl(argv[3], DIMS, in_dims); int dim = atoi(argv[1]); int rep = atoi(argv[2]); assert(dim < DIMS); assert(rep >= 0); assert(1 == in_dims[dim]); md_copy_dims(DIMS, out_dims, in_dims); out_dims[dim] = rep; complex float* out_data = create_cfl(argv[4], DIMS, out_dims); long in_strs[DIMS]; long out_strs[DIMS]; md_calc_strides(DIMS, in_strs, in_dims, CFL_SIZE); md_calc_strides(DIMS, out_strs, out_dims, CFL_SIZE); md_copy2(DIMS, out_dims, out_strs, out_data, in_strs, in_data, CFL_SIZE); unmap_cfl(DIMS, out_dims, out_data); unmap_cfl(DIMS, in_dims, in_data); exit(0); }
/** * Create a linear operator (without strides) * * @param N number of dimensions * @param odims dimensions of output (codomain) * @param idims dimensions of input (domain) * @param data data for applying the operator * @param forward function for applying the forward operation, A * @param adjoint function for applying the adjoint operation, A^H * @param normal function for applying the normal equations operation, A^H A * @param norm_inv function for applying the pseudo-inverse operation, (A^H A + mu I)^-1 * @param del function for freeing the data */ struct linop_s* linop_create(unsigned int ON, const long odims[ON], unsigned int IN, const long idims[IN], void* data, lop_fun_t forward, lop_fun_t adjoint, lop_fun_t normal, lop_p_fun_t norm_inv, del_fun_t del) { long ostrs[ON]; long istrs[IN]; md_calc_strides(ON, ostrs, odims, CFL_SIZE); md_calc_strides(IN, istrs, idims, CFL_SIZE); return linop_create2(ON, odims, ostrs, IN, idims, istrs, data, forward, adjoint, normal, norm_inv, del); }
void fwtN(unsigned int N, unsigned int flags, const long shifts[N], const long dims[N], const long ostr[2 * N], complex float* out, const long istr[N], const complex float* in, const long flen, const float filter[2][2][flen]) { long odims[2 * N]; wavelet_dims(N, flags, odims, dims, flen); assert(md_calc_size(2 * N, odims) >= md_calc_size(N, dims)); // FIXME one of these is unnecessary if we use the output complex float* tmpA = md_alloc_sameplace(2 * N, odims, CFL_SIZE, out); complex float* tmpB = md_alloc_sameplace(2 * N, odims, CFL_SIZE, out); long tidims[2 * N]; md_copy_dims(N, tidims, dims); md_singleton_dims(N, tidims + N); long tistrs[2 * N]; md_calc_strides(2 * N, tistrs, tidims, CFL_SIZE); long todims[2 * N]; md_copy_dims(2 * N, todims, tidims); long tostrs[2 * N]; // maybe we should push the randshift into lower levels //md_copy2(N, dims, tistrs, tmpA, istr, in, CFL_SIZE); md_circ_shift2(N, dims, shifts, tistrs, tmpA, istr, in, CFL_SIZE); for (unsigned int i = 0; i < N; i++) { if (MD_IS_SET(flags, i)) { todims[0 + i] = odims[0 + i]; todims[N + i] = odims[N + i]; md_calc_strides(2 * N, tostrs, todims, CFL_SIZE); fwt1(2 * N, i, tidims, tostrs, tmpB, (void*)tmpB + tostrs[N + i], tistrs, tmpA, flen, filter); md_copy_dims(2 * N, tidims, todims); md_copy_dims(2 * N, tistrs, tostrs); complex float* swap = tmpA; tmpA = tmpB; tmpB = swap; } } md_copy2(2 * N, todims, ostr, out, tostrs, tmpA, CFL_SIZE); md_free(tmpA); md_free(tmpB); }
struct linop_s* sampling_create(const long dims[DIMS], const long pat_dims[DIMS], const complex float* pattern) { struct sampling_data_s* data = xmalloc(sizeof(struct sampling_data_s)); md_select_dims(DIMS, ~MAPS_FLAG, data->dims, dims); // dimensions of kspace md_calc_strides(DIMS, data->strs, data->dims, CFL_SIZE); md_calc_strides(DIMS, data->pat_strs, pat_dims, CFL_SIZE); data->pattern = pattern; return linop_create(DIMS, data->dims, DIMS, data->dims, data, sampling_apply, sampling_apply, sampling_apply, NULL, sampling_free); }
void iwtN(unsigned int N, unsigned int flags, const long shifts[N], const long dims[N], const long ostr[N], complex float* out, const long istr[2 * N], const complex float* in, const long flen, const float filter[2][2][flen]) { long idims[2 * N]; wavelet_dims(N, flags, idims, dims, flen); assert(md_calc_size(2 * N, idims) >= md_calc_size(N, dims)); complex float* tmpA = md_alloc_sameplace(2 * N, idims, CFL_SIZE, out); complex float* tmpB = md_alloc_sameplace(2 * N, idims, CFL_SIZE, out); long tidims[2 * N]; md_copy_dims(2 * N, tidims, idims); long tistrs[2 * N]; md_calc_strides(2 * N, tistrs, tidims, CFL_SIZE); long todims[2 * N]; md_copy_dims(2 * N, todims, tidims); long tostrs[2 * N]; long ishifts[N]; for (unsigned int i = 0; i < N; i++) ishifts[i] = -shifts[i]; md_copy2(2 * N, tidims, tistrs, tmpA, istr, in, CFL_SIZE); for (int i = N - 1; i >= 0; i--) { // run backwards to maintain contigous blocks if (MD_IS_SET(flags, i)) { todims[0 + i] = dims[0 + i]; todims[N + i] = 1; md_calc_strides(2 * N, tostrs, todims, CFL_SIZE); iwt1(2 * N, i, todims, tostrs, tmpB, tistrs, tmpA, (void*)tmpA + tistrs[N + i], flen, filter); md_copy_dims(2 * N, tidims, todims); md_copy_dims(2 * N, tistrs, tostrs); complex float* swap = tmpA; tmpA = tmpB; tmpB = swap; } } //md_copy2(N, dims, ostr, out, tostrs, tmpA, CFL_SIZE); md_circ_shift2(N, dims, ishifts, ostr, out, tostrs, tmpA, CFL_SIZE); md_free(tmpA); md_free(tmpB); }
static complex float* compute_psf2(unsigned int N, const long psf_dims[N + 3], const long trj_dims[N + 3], const complex float* traj, const complex float* weights) { unsigned int ND = N + 3; long img_dims[ND]; long img_strs[ND]; md_select_dims(ND, ~MD_BIT(N + 0), img_dims, psf_dims); md_calc_strides(ND, img_strs, img_dims, CFL_SIZE); // PSF 2x size long img2_dims[ND]; long img2_strs[ND]; md_copy_dims(ND, img2_dims, img_dims); for (int i = 0; i < 3; i++) img2_dims[i] = (1 == img_dims[i]) ? 1 : (2 * img_dims[i]); md_calc_strides(ND, img2_strs, img2_dims, CFL_SIZE); complex float* traj2 = md_alloc(ND, trj_dims, CFL_SIZE); md_zsmul(ND, trj_dims, traj2, traj, 2.); complex float* psft = compute_psf(ND, img2_dims, trj_dims, traj2, weights); md_free(traj2); fftuc(ND, img2_dims, FFT_FLAGS, psft, psft); // reformat long sub2_strs[ND]; md_copy_strides(ND, sub2_strs, img2_strs); for(int i = 0; i < 3; i++) sub2_strs[i] *= 2;; complex float* psf = md_alloc(ND, psf_dims, CFL_SIZE); long factors[N]; for (unsigned int i = 0; i < N; i++) factors[i] = ((img_dims[i] > 1) && (i < 3)) ? 2 : 1; md_decompose(N + 0, factors, psf_dims, psf, img2_dims, psft, CFL_SIZE); md_free(psft); return psf; }
static double bench_generic_add(long dims[DIMS], unsigned int flags, bool forloop) { long dimsX[DIMS]; long dimsY[DIMS]; long dimsC[DIMS]; md_select_dims(DIMS, flags, dimsX, dims); md_select_dims(DIMS, ~flags, dimsC, dims); md_select_dims(DIMS, ~0u, dimsY, dims); long strsX[DIMS]; long strsY[DIMS]; md_calc_strides(DIMS, strsX, dimsX, CFL_SIZE); md_calc_strides(DIMS, strsY, dimsY, CFL_SIZE); complex float* x = md_alloc(DIMS, dimsX, CFL_SIZE); complex float* y = md_alloc(DIMS, dimsY, CFL_SIZE); md_gaussian_rand(DIMS, dimsX, x); md_gaussian_rand(DIMS, dimsY, y); long L = md_calc_size(DIMS, dimsC); long T = md_calc_size(DIMS, dimsX); double tic = timestamp(); if (forloop) { for (long i = 0; i < L; i++) { for (long j = 0; j < T; j++) y[i + j * L] += x[j]; } } else { md_zaxpy2(DIMS, dims, strsY, y, 1., strsX, x); } double toc = timestamp(); md_free(x); md_free(y); return toc - tic; }
int main_mip(int argc, char* argv[argc]) { bool mIP = false; const struct opt_s opts[] = { OPT_SET('m', &mIP, "minimum" ), }; cmdline(&argc, argv, 3, 3, usage_str, help_str, ARRAY_SIZE(opts), opts); unsigned int flags = atoi(argv[1]); long idims[DIMS]; complex float* in = load_cfl(argv[2], DIMS, idims); long odims[DIMS]; md_select_dims(DIMS, ~flags, odims, idims); complex float* out = create_cfl(argv[3], DIMS, odims); complex float* tmp = md_alloc(DIMS, idims, CFL_SIZE); md_zabs(DIMS, idims, tmp, in); long istr[DIMS]; long ostr[DIMS]; md_calc_strides(DIMS, istr, idims, CFL_SIZE); md_calc_strides(DIMS, ostr, odims, CFL_SIZE); md_clear(DIMS, odims, out, CFL_SIZE); md_max2(DIMS, idims, ostr, (float*)out, ostr, (const float*)out, istr, (const float*)tmp); if (mIP) { // need result of max in output md_min2(DIMS, idims, ostr, (float*)out, ostr, (const float*)out, istr, (const float*)tmp); } md_free(tmp); unmap_cfl(DIMS, idims, in); unmap_cfl(DIMS, odims, out); exit(0); }
struct prox_dfwavelet_data* prepare_prox_dfwavelet_data(const long im_dims[DIMS], const long min_size[3], const complex float res[3], unsigned int flow_dim, float lambda, bool use_gpu) { // get dimension PTR_ALLOC(struct prox_dfwavelet_data, data); md_copy_dims(DIMS, data->im_dims, im_dims); md_select_dims(DIMS, FFT_FLAGS, data->tim_dims, im_dims); md_calc_strides(DIMS, data->im_strs, im_dims, CFL_SIZE); // initialize temp #ifdef USE_CUDA if (use_gpu) { data->vx = md_alloc_gpu(DIMS, data->tim_dims, CFL_SIZE); data->vy = md_alloc_gpu(DIMS, data->tim_dims, CFL_SIZE); data->vz = md_alloc_gpu(DIMS, data->tim_dims, CFL_SIZE); } else #endif { data->vx = md_alloc(DIMS, data->tim_dims, CFL_SIZE); data->vy = md_alloc(DIMS, data->tim_dims, CFL_SIZE); data->vz = md_alloc(DIMS, data->tim_dims, CFL_SIZE); } data->flow_dim = flow_dim; data->slice_flag = ~FFT_FLAGS; data->lambda = lambda; data->plan = prepare_dfwavelet_plan(3, data->tim_dims, (long*) min_size, (complex float*) res, use_gpu); return data; }
void fftscale(unsigned int N, const long dims[N], unsigned long flags, complex float* dst, const complex float* src) { long strs[N]; md_calc_strides(N, strs, dims, CFL_SIZE); fftscale2(N, dims, flags, strs, dst, strs, src); }
const struct operator_s* fft_create(unsigned int D, const long dimensions[D], unsigned long flags, complex float* dst, const complex float* src, bool backwards) { long strides[D]; md_calc_strides(D, strides, dimensions, CFL_SIZE); return fft_create2(D, dimensions, flags, strides, dst, strides, src, backwards); }
/** * Thresholding operator for l0-norm: f(x) = || x ||_0 <= k, as used in NIHT algorithm. * y = HT(x, k) (hard thresholding, ie keeping the k largest elements). * * @param D number of dimensions * @param dim dimensions of x * @param k threshold parameter (non-zero elements to keep) * @param flags bitmask for joint thresholding */ const struct operator_p_s* prox_niht_thresh_create(unsigned int D, const long dim[D], const unsigned int k, const unsigned long flags) { PTR_ALLOC(struct thresh_s, data); SET_TYPEID(thresh_s, data); data->lambda = 0.; data->k = k; data->D = D; data->flags = flags; data->unitary_op = NULL; PTR_ALLOC(long[D], ndim); md_copy_dims(D, *ndim, dim); data->dim = *PTR_PASS(ndim); // norm dimensions are the flagged input dimensions PTR_ALLOC(long[D], norm_dim); md_select_dims(D, ~flags, *norm_dim, data->dim); data->norm_dim = *PTR_PASS(norm_dim); PTR_ALLOC(long[D], nstr); md_calc_strides(D, *nstr, data->dim, CFL_SIZE); data->str = *PTR_PASS(nstr); return operator_p_create(D, dim, D, dim, CAST_UP(PTR_PASS(data)), hardthresh_apply, thresh_del); }
/** * Proximal operator for l1-norm with unitary transform: f(x) = lambda || T x ||_1 * * @param D number of dimensions * @param dim dimensions of x * @param lambda threshold parameter * @param unitary_op unitary linear operator * @param flags bitmask for joint soft-thresholding */ extern const struct operator_p_s* prox_unithresh_create(unsigned int D, const struct linop_s* unitary_op, const float lambda, const unsigned long flags) { PTR_ALLOC(struct thresh_s, data); SET_TYPEID(thresh_s, data); data->lambda = lambda; data->D = D; data->flags = flags; data->unitary_op = unitary_op; const long* dims = linop_domain(unitary_op)->dims; PTR_ALLOC(long[D], ndim); md_copy_dims(D, *ndim, dims); data->dim = *PTR_PASS(ndim); PTR_ALLOC(long[D], nstr); md_calc_strides(D, *nstr, data->dim, CFL_SIZE); data->str = *PTR_PASS(nstr); // norm dimensions are the flagged transform dimensions // FIXME should use linop_codomain(unitary_op)->N PTR_ALLOC(long[D], norm_dim); md_select_dims(D, ~flags, *norm_dim, linop_codomain(unitary_op)->dims); data->norm_dim = *PTR_PASS(norm_dim); return operator_p_create(D, dims, D, dims, CAST_UP(PTR_PASS(data)), unisoftthresh_apply, thresh_del); }
void fwt(unsigned int N, unsigned int flags, const long shifts[N], const long dims[N], complex float* out, const long istr[N], const complex float* in, const long minsize[N], long flen, const float filter[2][2][flen]) { if (0 == flags) { if (out != in) md_copy2(N, dims, istr, out, istr, in, CFL_SIZE); return; } unsigned long coeffs = wavelet_coeffs(N, flags, dims, minsize, flen); long wdims[2 * N]; wavelet_dims(N, flags, wdims, dims, flen); long ostr[2 * N]; md_calc_strides(2 * N, ostr, wdims, CFL_SIZE); long offset = coeffs - md_calc_size(2 * N, wdims); debug_printf(DP_DEBUG4, "%d %ld %ld\n", flags, coeffs, offset); long shifts0[N]; for (unsigned int i = 0; i < N; i++) shifts0[i] = 0; fwtN(N, flags, shifts, dims, ostr, out + offset, istr, in, flen, filter); fwt(N, wavelet_filter_flags(N, flags, wdims, minsize), shifts0, wdims, out, ostr, out + offset, minsize, flen, filter); }
int main_reshape(int argc, char* argv[]) { cmdline(&argc, argv, 3, 100, usage_str, help_str, 0, NULL); num_init(); unsigned int flags = atoi(argv[1]); unsigned int n = bitcount(flags); assert((int)n + 3 == argc - 1); long in_dims[DIMS]; long in_strs[DIMS]; long out_dims[DIMS]; long out_strs[DIMS]; complex float* in_data = load_cfl(argv[n + 2], DIMS, in_dims); md_calc_strides(DIMS, in_strs, in_dims, CFL_SIZE); md_copy_dims(DIMS, out_dims, in_dims); unsigned int j = 0; for (unsigned int i = 0; i < DIMS; i++) if (MD_IS_SET(flags, i)) out_dims[i] = atoi(argv[j++ + 2]); assert(j == n); assert(md_calc_size(DIMS, in_dims) == md_calc_size(DIMS, out_dims)); md_calc_strides(DIMS, out_strs, out_dims, CFL_SIZE); for (unsigned int i = 0; i < DIMS; i++) if (!(MD_IS_SET(flags, i) || (in_strs[i] == out_strs[i]))) error("Dimensions are not consistent at index %d.\n"); complex float* out_data = create_cfl(argv[n + 3], DIMS, out_dims); md_copy(DIMS, in_dims, out_data, in_data, CFL_SIZE); unmap_cfl(DIMS, in_dims, in_data); unmap_cfl(DIMS, out_dims, out_data); exit(0); }
static struct linop_s* linop_fft_create_priv(int N, const long dims[N], unsigned int flags, bool forward, bool center) { const struct operator_s* plan = fft_measure_create(N, dims, flags, true, false); const struct operator_s* iplan = fft_measure_create(N, dims, flags, true, true); PTR_ALLOC(struct fft_linop_s, data); SET_TYPEID(fft_linop_s, data); data->frw = plan; data->adj = iplan; data->N = N; data->center = center; data->dims = *TYPE_ALLOC(long[N]); md_copy_dims(N, data->dims, dims); data->strs = *TYPE_ALLOC(long[N]); md_calc_strides(N, data->strs, data->dims, CFL_SIZE); long fft_dims[N]; md_select_dims(N, flags, fft_dims, dims); data->nscale = (float)md_calc_size(N, fft_dims); lop_fun_t apply = forward ? fft_linop_apply : fft_linop_adjoint; lop_fun_t adjoint = forward ? fft_linop_adjoint : fft_linop_apply; struct linop_s* lop = linop_create(N, dims, N, dims, CAST_UP(PTR_PASS(data)), apply, adjoint, fft_linop_normal, NULL, fft_linop_free); if (center) { // FIXME: should only allocate flagged dims complex float* fftmod_mat = md_alloc(N, dims, CFL_SIZE); complex float* fftmodk_mat = md_alloc(N, dims, CFL_SIZE); // we need fftmodk only because we want to apply scaling only once complex float one[1] = { 1. }; md_fill(N, dims, fftmod_mat, one, CFL_SIZE); fftmod(N, dims, flags, fftmodk_mat, fftmod_mat); fftscale(N, dims, flags, fftmod_mat, fftmodk_mat); struct linop_s* mod = linop_cdiag_create(N, dims, ~0u, fftmod_mat); struct linop_s* modk = linop_cdiag_create(N, dims, ~0u, fftmodk_mat); struct linop_s* tmp = linop_chain(mod, lop); tmp = linop_chain(tmp, modk); linop_free(lop); linop_free(mod); linop_free(modk); lop = tmp; } return lop; }
void overlapandsave(int N, const long dims[N], const long blk[N], complex float* dst, complex float* src1, const long dim2[N], complex float* src2) { // [------++++ // [------ long ndims[2 * N]; long L[2 * N]; long ndim2[2 * N]; long ndim3[2 * N]; for (int i = 0; i < N; i++) { assert(0 == dims[i] % blk[i]); assert(dim2[i] <= blk[i]); ndims[i * 2 + 1] = dims[i] / blk[i]; ndims[i * 2 + 0] = blk[i]; L[i * 2 + 1] = dims[i] / blk[i]; L[i * 2 + 0] = blk[i] + dim2[i] - 1; ndim2[i * 2 + 1] = 1; ndim2[i * 2 + 0] = dim2[i]; ndim3[i * 2 + 1] = dims[i] / blk[i] - 0; ndim3[i * 2 + 0] = blk[i]; } long T = md_calc_size(2 * N, L); complex float* tmp = xmalloc(T * 8); long str1[2 * N]; long str2[2 * N]; long str3[2 * N]; md_calc_strides(2 * N, str1, ndims, 8); md_calc_strides(2 * N, str2, L, 8); md_calc_strides(2 * N, str3, ndim3, 8); md_clear(2 * N, L, tmp, 8); md_copy2(2 * N, ndim3, str2, tmp, str1, src1, 8); conv(2 * N, ~0, CONV_VALID, CONV_CAUSAL, ndims, dst, L, tmp, ndim2, src2); free(tmp); }
void overlapandadd(int N, const long dims[N], const long blk[N], complex float* dst, complex float* src1, const long dim2[N], complex float* src2) { long ndims[2 * N]; long L[2 * N]; long ndim2[2 * N]; long ndim3[2 * N]; for (int i = 0; i < N; i++) { assert(0 == dims[i] % blk[i]); assert(dim2[i] <= blk[i]); ndims[i * 2 + 1] = dims[i] / blk[i]; ndims[i * 2 + 0] = blk[i]; L[i * 2 + 1] = dims[i] / blk[i]; L[i * 2 + 0] = blk[i] + dim2[i] - 1; ndim2[i * 2 + 1] = 1; ndim2[i * 2 + 0] = dim2[i]; ndim3[i * 2 + 1] = dims[i] / blk[i] + 1; ndim3[i * 2 + 0] = blk[i]; } long T = md_calc_size(2 * N, L); complex float* tmp = xmalloc(T * 8); // conv_causal_extend(2 * N, L, tmp, ndims, src1, ndim2, src2); conv(2 * N, ~0, CONV_EXTENDED, CONV_CAUSAL, L, tmp, ndims, src1, ndim2, src2); // [------++++|||||||| //long str1[2 * N]; long str2[2 * N]; long str3[2 * N]; //md_calc_strides(2 * N, str1, ndims, 8); md_calc_strides(2 * N, str2, L, 8); md_calc_strides(2 * N, str3, ndim3, 8); md_clear(2 * N, ndim3, dst, CFL_SIZE); md_zadd2(2 * N, L, str3, dst, str3, dst, str2, tmp); free(tmp); }
struct linop_s* linop_sampling_create(const long dims[DIMS], const long pat_dims[DIMS], const complex float* pattern) { PTR_ALLOC(struct sampling_data_s, data); SET_TYPEID(sampling_data_s, data); md_copy_dims(DIMS, data->pat_dims, pat_dims); md_select_dims(DIMS, ~MAPS_FLAG, data->dims, dims); // dimensions of kspace md_calc_strides(DIMS, data->strs, data->dims, CFL_SIZE); md_calc_strides(DIMS, data->pat_strs, data->pat_dims, CFL_SIZE); data->pattern = (complex float*)pattern; #ifdef USE_CUDA data->gpu_pattern = NULL; #endif const long* dims2 = data->dims; return linop_create(DIMS, dims2, DIMS, dims2, CAST_UP(PTR_PASS(data)), sampling_apply, sampling_apply, sampling_apply, NULL, sampling_free); }
/* * Nuclear norm calculation for arbitrary block sizes */ float lrnucnorm(const struct operator_p_s* op, const complex float* src) { struct lrthresh_data_s* data = (struct lrthresh_data_s*)operator_p_get_data(op); long strs1[DIMS]; md_calc_strides(DIMS, strs1, data->dims_decom, 1); float nnorm = 0.; for (int l = 0; l < data->levels; l++) { const complex float* srcl = src + l * strs1[LEVEL_DIM]; // Initialize long blkdims[DIMS]; long blksize = 1; for (unsigned int i = 0; i < DIMS; i++) { blkdims[i] = data->blkdims[l][i]; blksize *= blkdims[i]; } // Special case if blocksize is 1 if (blksize == 1) { for (long j = 0; j < md_calc_size(DIMS, data->dims); j++) nnorm += 2 * cabsf(srcl[j]); continue; } // Initialize data struct svthresh_blockproc_data* svdata = svthresh_blockproc_create(data->mflags, 0., 0); // Initialize tmp complex float* tmp; #ifdef USE_CUDA tmp = (data->use_gpu ? md_alloc_gpu : md_alloc)(DIMS, data->dims, CFL_SIZE); #else tmp = md_alloc(DIMS, data->dims, CFL_SIZE); #endif // Copy to tmp //debug_print_dims(DP_DEBUG1, DIMS, data->dims); md_copy(DIMS, data->dims, tmp, srcl, CFL_SIZE); // Block SVD Threshold nnorm = blockproc(DIMS, data->dims, blkdims, (void*)svdata, nucnorm_blockproc, tmp, tmp); // Free tmp free(svdata); md_free(tmp); } return nnorm; }
const struct linop_s* linop_fmac_create(unsigned int N, const long dims[N], unsigned int oflags, unsigned int iflags, unsigned int tflags, const complex float* tensor) { PTR_ALLOC(struct fmac_data, data); SET_TYPEID(fmac_data, data); data->N = N; data->dims = *TYPE_ALLOC(long[N]); md_copy_dims(N, data->dims, dims); data->idims = *TYPE_ALLOC(long[N]); data->istrs = *TYPE_ALLOC(long[N]); md_select_dims(N, ~iflags, data->idims, dims); md_calc_strides(N, data->istrs, data->idims, CFL_SIZE); data->odims = *TYPE_ALLOC(long[N]); data->ostrs = *TYPE_ALLOC(long[N]); md_select_dims(N, ~oflags, data->odims, dims); md_calc_strides(N, data->ostrs, data->odims, CFL_SIZE); data->tstrs = *TYPE_ALLOC(long[N]); data->tdims = *TYPE_ALLOC(long[N]); md_select_dims(N, ~tflags, data->tdims, dims); md_calc_strides(N, data->tstrs, data->tdims, CFL_SIZE); data->tensor = tensor; #ifdef USE_CUDA data->gpu_tensor = NULL; #endif long odims[N]; md_copy_dims(N, odims, data->odims); long idims[N]; md_copy_dims(N, idims, data->idims); return linop_create(N, odims, N, idims, CAST_UP(PTR_PASS(data)), fmac_apply, fmac_adjoint, NULL, NULL, fmac_free_data); }
void data_consistency(const long dims[DIMS], complex float* dst, const complex float* pattern, const complex float* kspace1, const complex float* kspace2) { assert(1 == dims[MAPS_DIM]); long strs[DIMS]; long dims1[DIMS]; long strs1[DIMS]; md_select_dims(DIMS, ~COIL_FLAG, dims1, dims); md_calc_strides(DIMS, strs1, dims1, CFL_SIZE); md_calc_strides(DIMS, strs, dims, CFL_SIZE); complex float* tmp = md_alloc_sameplace(DIMS, dims, CFL_SIZE, dst); md_zmul2(DIMS, dims, strs, tmp, strs, kspace2, strs1, pattern); md_zsub(DIMS, dims, tmp, kspace2, tmp); md_zfmac2(DIMS, dims, strs, tmp, strs, kspace1, strs1, pattern); md_copy(DIMS, dims, dst, tmp, CFL_SIZE); md_free(tmp); }
/** * Singular Value Thresholding * * @param M - matrix column size * @param N - matrix row size * @param lambda - regularization parameter * @param A - input/output matrix */ float svthresh(long M, long N, float lambda, complex float* dst, const complex float* src) //FIXME: destroys input { long minMN = MIN(M,N); long dimsU[3] = {M,minMN,1}; long dimsVT[3] = {minMN,N,1}; long dimsS[3] = {minMN,1,1}; // long dimsAA[3] = {minMN, minMN,1}; long strsVT[3]; long strsS[3]; md_calc_strides(3, strsVT, dimsVT, CFL_SIZE); md_calc_strides(3, strsS, dimsS, FL_SIZE); complex float* U = md_alloc_sameplace(3, dimsU, CFL_SIZE, src); complex float* VT = md_alloc_sameplace(3, dimsVT, CFL_SIZE, src ); float* S = md_alloc_sameplace(3, dimsS, FL_SIZE, src); // complex float* AA = md_alloc_sameplace(3, dimsAA, CFL_SIZE, src ); // lapack_normal_multiply( M, N, (M > N), (complex float (*) [])AA, (const complex float (*) [])src ); // SVD lapack_svd_econ(M, N, (complex float (*) []) U, (complex float (*) []) VT, S, (complex float (*) [N])src); // Thresh S md_softthresh(3, dimsS, lambda, 0, S, S); // VT = S * VT md_mul2( 3, dimsVT, strsVT, (float*) VT, strsVT, (float*) VT, strsS, S ); md_mul2( 3, dimsVT, strsVT, ((float*) VT)+1, strsVT, ((float*) VT)+1, strsS, S ); // dst = U * VT blas_matrix_multiply( M, N, minMN, (complex float (*) [N])dst, (const complex float (*) [minMN])U, (const complex float (*) [N])VT ); md_free(U); md_free(VT); md_free(S); return 0; }
static struct sum_data* sum_create_data( const long imgd_dims[DIMS], bool use_gpu ) { PTR_ALLOC(struct sum_data, data); // decom dimensions md_copy_dims(DIMS, data->imgd_dims, imgd_dims); md_calc_strides(DIMS, data->imgd_strs, imgd_dims, CFL_SIZE); // image dimensions data->levels = imgd_dims[LEVEL_DIM]; md_select_dims(DIMS, ~LEVEL_FLAG, data->img_dims, imgd_dims); md_calc_strides(DIMS, data->img_strs, data->img_dims, CFL_SIZE); data->tmp = NULL; data->use_gpu = use_gpu; return data; }