// Declarations of classes that implements useful // utility functions for the eddy current project. // They are collections of statically declared // functions that have been collected into classes // to make it explicit where they come from. There // will never be any instances of theses classes. // // EddyUtils.cpp // // Jesper Andersson, FMRIB Image Analysis Group // // Copyright (C) 2011 University of Oxford // #include #include #include #include #include #include "armawrap/newmat.h" #ifndef EXPOSE_TREACHEROUS #define EXPOSE_TREACHEROUS // To allow us to use .set_sform etc #endif #include "newimage/newimageall.h" #include "miscmaths/miscmaths.h" #include "cprob/libprob.h" #include "warpfns/warpfns.h" #include "EddyHelperClasses.h" #include "EddyUtils.h" using namespace std; using namespace CPROB; using namespace EDDY; // Remove this line when we move NVCC to C++17 double EddyUtils::b_range = 100; //@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@ // // Class EddyUtils // // Helper Class used to perform various useful tasks for // the eddy current correction project. // //@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@ bool EddyUtils::get_groups(// Input const std::vector& dpv, // Output std::vector >& grps, std::vector& grpi, std::vector& grpb) EddyTry { std::vector grp_templates; // First pass to sort out how many different b-values/shells there are grp_templates.push_back(0); for (unsigned int i=1; i grp_n(grp_templates.size(),1); for (unsigned int j=0; j((double(dpv.size() - grp_n[0]) / double(grpb.size() - 1)) + 0.5); is_shelled &= bool(*std::max_element(grp_n.begin()+1,grp_n.end()) < 2 * scans_per_shell); // Don't trust too many scans in one shell is_shelled &= bool(3 * *std::min_element(grp_n.begin()+1,grp_n.end()) > scans_per_shell); // Don't trust too few scans in one shell } else { // If all scans are dwis is_shelled = grpb.size() < 6; // Don't trust more than 5 shells unsigned int scans_per_shell = static_cast((double(dpv.size()) / double(grpb.size())) + 0.5); is_shelled &= bool(*std::max_element(grp_n.begin(),grp_n.end()) < 2 * scans_per_shell); // Don't trust too many scans in one shell is_shelled &= bool(3 * *std::min_element(grp_n.begin(),grp_n.end()) > scans_per_shell); // Don't trust too few scans in one shell } if (!is_shelled) return(false); // Final sanity check unsigned int nscan = grps[0].size(); for (unsigned int i=1; i& dpv) EddyTry { std::vector > grps; std::vector grpi; std::vector grpb; return(get_groups(dpv,grps,grpi,grpb)); } EddyCatch bool EddyUtils::IsMultiShell(const std::vector& dpv) EddyTry { std::vector > grps; std::vector grpi; std::vector grpb; bool is_shelled = get_groups(dpv,grps,grpi,grpb); return(is_shelled && grpb.size() > 1); } EddyCatch bool EddyUtils::GetGroups(// Input const std::vector& dpv, // Output std::vector& grpi, std::vector& grpb) EddyTry { std::vector > grps; return(get_groups(dpv,grps,grpi,grpb)); } EddyCatch bool EddyUtils::GetGroups(// Input const std::vector& dpv, // Output std::vector >& grps, std::vector& grpb) EddyTry { std::vector grpi; return(get_groups(dpv,grps,grpi,grpb)); } EddyCatch bool EddyUtils::GetGroups(// Input const std::vector& dpv, // Output std::vector >& grps, std::vector& grpi, std::vector& grpb) EddyTry { return(get_groups(dpv,grps,grpi,grpb)); } EddyCatch std::string EddyUtils::ModelString(ECModelType mt) EddyTry { switch (mt) { case ECModelType::NoEC: return("No EC model"); case ECModelType::Linear: return("Linear EC model"); case ECModelType::Quadratic: return("Quadratic EC model"); case ECModelType::Cubic: return("Cubic EC model"); default: return("Unknown EC model"); } } EddyCatch std::vector EddyUtils::GetIndiciesOfDWIs(const std::vector& dpars) EddyTry { std::vector indicies; for (unsigned int i=0; i EddyUtils::GetDWIDiffParas(const std::vector& dpars) EddyTry { std::vector indx = EddyUtils::GetIndiciesOfDWIs(dpars); std::vector dwi_dpars; for (unsigned int i=0; i 1e-6) return(false); if (!EddyUtils::AreInSameShell(s1.GetDiffPara(),s2.GetDiffPara())) return(false); if (IsDiffusionWeighted(s1.GetDiffPara()) && !HaveSameDirection(s1.GetDiffPara(),s2.GetDiffPara())) return(false); return(true); } EddyCatch std::vector EddyUtils::GetSliceWiseForwardMovementMatrices(const EDDY::ECScan& scan) EddyTry { std::vector R(scan.GetIma().zsize()); for (unsigned int tp=0; tp slices = scan.GetMBG().SlicesAtTimePoint(tp); for (unsigned int i=0; i EddyUtils::GetSliceWiseInverseMovementMatrices(const EDDY::ECScan& scan) EddyTry { std::vector R(scan.GetIma().zsize()); for (unsigned int tp=0; tp slices = scan.GetMBG().SlicesAtTimePoint(tp); for (unsigned int i=0; i& dwivols, const std::string& fname, const std::vector& dpars) EddyTry { std::vector indx = EddyUtils::GetIndiciesOfDWIs(dpars); NEWIMAGE::volume tmp; read_volumeROI(tmp,fname,0,0,0,0,-1,-1,-1,0); dwivols.reinitialize(tmp.xsize(),tmp.ysize(),tmp.zsize(),indx.size()); for (unsigned int i=0; i EddyUtils::ConvertMaskToFloat(const NEWIMAGE::volume& charmask) EddyTry { NEWIMAGE::volume floatmask(charmask.xsize(),charmask.ysize(),charmask.zsize()); NEWIMAGE::copybasicproperties(charmask,floatmask); for (int k=0; k(charmask(i,j,k)); } } } return(floatmask); } EddyCatch // Rewritten for new newimage NEWIMAGE::volume EddyUtils::Smooth(const NEWIMAGE::volume& ima, float fwhm, const NEWIMAGE::volume& mask) EddyTry { if (mask.getextrapolationmethod() != NEWIMAGE::zeropad) throw EddyException("EddyUtils::Smooth: mask must use zeropad for extrapolation"); float sx = (fwhm/std::sqrt(8.0*std::log(2.0)))/ima.xdim(); float sy = (fwhm/std::sqrt(8.0*std::log(2.0)))/ima.ydim(); float sz = (fwhm/std::sqrt(8.0*std::log(2.0)))/ima.zdim(); int nx=((int) (sx-0.001))*2 + 3; int ny=((int) (sy-0.001))*2 + 3; int nz=((int) (sz-0.001))*2 + 3; NEWMAT::ColumnVector krnlx = NEWIMAGE::gaussian_kernel1D(sx,nx); NEWMAT::ColumnVector krnly = NEWIMAGE::gaussian_kernel1D(sy,ny); NEWMAT::ColumnVector krnlz = NEWIMAGE::gaussian_kernel1D(sz,nz); NEWIMAGE::volume4D ovol = ima; // volume4D just an alias for volume for (int i=0; i EddyUtils::Smooth4D(const NEWIMAGE::volume4D& ima, float fwhm, const NEWIMAGE::volume& mask) { NEWIMAGE::volume4D ovol = ima; for (int i=0; i EddyUtils::Negate(const std::vector& vec) EddyTry { std::vector rvec = vec; std::for_each(rvec.begin(),rvec.end(),[](arma::colvec& v){ v = -v; }); return(rvec); } EddyCatch NEWIMAGE::volume EddyUtils::MakeNoiseIma(const NEWIMAGE::volume& ima, float mu, float stdev) EddyTry { NEWIMAGE::volume nima = ima; double rnd; for (int k=0; k(ndtri(rnd-1)); } } } return(nima); } EddyCatch DiffStats EddyUtils::GetSliceWiseStats(// Input const EddyCommandLineOptions& clo, ScanType st, const NEWIMAGE::volume& pred, // Prediction in model space std::shared_ptr > susc, // Susceptibility induced off-resonance field const NEWIMAGE::volume& pmask, // "Data valid" mask in model space const ECScanManager& sm, unsigned int s, // Scan index unsigned int iter, // Iteration unsigned int dl) EddyTry // Debug level { // Transform prediction into observation space NEWIMAGE::volume mask = pred; mask = 0.0; NEWIMAGE::volume jac; NEWIMAGE::volume pios = EddyUtils::transform_model_to_scan_space(pred,sm.Scan(s,st),susc,true,mask,&jac,nullptr); // Transform binary mask into observation space mask = 0.0; NEWIMAGE::volume bios = EddyUtils::transform_model_to_scan_space(pmask*sm.Mask(),sm.Scan(s,st),susc,false,mask,nullptr,nullptr); bios.binarise(0.99); // Value above (arbitrary) 0.99 implies valid voxels mask *= bios; // Volume and prediction mask falls within FOV jac.binarise(0.0,clo.OLUpperThresholdJacobianMask()); mask *= jac; // Calculate slice-wise stats from difference image DiffStats stats(sm.Scan(s,st).GetOriginalIma()-pios,mask); // Write debug info if requested if (dl > 2 && clo.DebugIndicies().IsAmongIndicies(sm.GetGlobalIndex(s,st))) { // Note that we use a local "debug level" instead of clo.DebugLevel() NEWIMAGE::volume out = sm.Scan(s,st).GetOriginalIma()-pios; out *= mask; std::string prefix("EDDY_DEBUG_"); if (st==ScanType::b0) prefix += "b0"; else if (st==ScanType::DWI) prefix += "DWI"; else if (st==ScanType::fMRI) prefix += "fMRI"; char fname[256]; sprintf(fname,"%s_repol_diff_%02d_%04d",prefix.c_str(),iter,sm.GetGlobalIndex(s,st)); NEWIMAGE::write_volume(out,fname); } return(stats); } EddyCatch std::vector EddyUtils::ScansPerThread(unsigned int nscans, // No. of scans/volumes unsigned int nthreads) EddyTry // No. of threads { double vpt = static_cast(nscans) / static_cast(nthreads); std::vector nvols(nthreads+1,0); for (unsigned int i=1; i(std::floor(i*vpt)); nvols[nthreads] = nscans; return(nvols); } EddyCatch double EddyUtils::param_update(// Input const NEWIMAGE::volume& pred, // Prediction in model space std::shared_ptr > susc, // Susceptibility off-resonance field std::shared_ptr > bias, // Recieve bias field const NEWIMAGE::volume& pmask, // "Data valid" mask in model space EDDY::ParametersType whichp, // Which parameters to update float fwhm, // FWHM for Gaussian smoothing bool very_verbose, // These input parameters are for debugging only unsigned int scindx, // Scan index unsigned int iter, // Iteration unsigned int level, // Determines how much gets written // Input/output EDDY::ECScan& scan) // Scan we want to register to pred EddyTry { static std::mutex cout_mutex; // Magic static // Transform prediction into observation space NEWIMAGE::volume mask = pred; mask.setextrapolationmethod(NEWIMAGE::zeropad); mask = 0.0; NEWIMAGE::volume jac = pred; jac = 1.0; NEWIMAGE::volume pios = EddyUtils::transform_model_to_scan_space(pred,scan,susc,true,mask,&jac,NULL); // NEWIMAGE::volume pios = EddyUtils::transform_model_to_scan_space(pred,scan,susc,true,mask,&jac,NULL,scindx,iter,level); // Transform binary mask into observation space NEWIMAGE::volume skrutt = pred; skrutt = 0.0; NEWIMAGE::volume mios = EddyUtils::transform_model_to_scan_space(pmask,scan,susc,false,skrutt,NULL,NULL); mios.binarise(0.99); // Value above (arbitrary) 0.99 implies valid voxels mask *= mios; // Volume and prediction mask falls within FOV // Get partial derivatives w.r.t. to requested category of parameters in prediction space NEWIMAGE::volume4D derivs = EddyUtils::get_partial_derivatives_in_scan_space(pred,scan,susc,whichp); if (fwhm) { mask.setextrapolationmethod(NEWIMAGE::zeropad); derivs = EddyUtils::Smooth(derivs,fwhm,mask); } // Calculate XtX where X is a matrix whose columns are the partial derivatives arma::mat XtX = EddyUtils::make_XtX(derivs,mask); // Calculate difference image between observed and predicted NEWIMAGE::volume dima = pios-scan.GetIma(); if (fwhm) { mask.setextrapolationmethod(NEWIMAGE::zeropad); dima = EddyUtils::Smooth(dima,fwhm,mask); } // Calculate Xty where y is the difference between observed and predicted. X as above. arma::colvec Xty = EddyUtils::make_Xty(derivs,dima,mask); // Get derivative and Hessian of regularisation (relevant only for slice-to-vol); arma::colvec lHb = scan.GetRegGrad(whichp); arma::mat H = scan.GetRegHess(whichp); // Calculate mean sum of squares from difference image and add regularisation double masksum = mask.sum(); double mss = (dima*mask).sumsquares() / masksum + scan.GetReg(whichp); // Very mild Tikhonov regularisation to select solution with smaller norm arma::mat regm = scan.GetLevenbergLambda() * arma::eye(XtX.n_rows,XtX.n_cols); // Calculate update to parameters arma::colvec update = -arma::solve(XtX/masksum+H+regm,Xty/masksum+lHb,arma::solve_opts::likely_sympd); // Calculate sims (scan in model space) if we need to write it as debug info NEWIMAGE::volume sims; if (level) sims = scan.GetUnwarpedIma(susc); // Update parameters for (unsigned int i=0; i(1,scindx)) + scan.GetReg(whichp); if (!std::isnan(new_mss) && new_mss < mss) { // Success scan.SetLevenbergLambda(0.1*scan.GetLevenbergLambda()); return(mss); } else { // mss didn't decrease. Increase regularisation until it does. while ((std::isnan(new_mss) || new_mss >= mss) && scan.GetLevenbergLambda() < 1.0e6) { if (very_verbose) { const std::lock_guard lock(cout_mutex); // Grab ownership of terminal output std::cout << "EddyUtils::param_update: updates rejected for volume #" << scindx << std::endl; std::cout << "EddyUtils::param_update: Levenberg-lambda = " << scan.GetLevenbergLambda(); std::cout << ", original mss = " << mss << ", after update mss = " << new_mss << std::endl; std::cout.flush(); } // Set back parameters to old values for (unsigned int i=0; i(1,scindx)) + scan.GetReg(whichp); } if (std::isnan(new_mss) || new_mss >= mss) { // It never decreased if (very_verbose) { const std::lock_guard lock(cout_mutex); // Grab ownership of terminal output std::cout << "EddyUtils::param_update: Final attempt, updates rejected for volume #" << scindx << std::endl; std::cout << "EddyUtils::param_update: Levenberg-lambda = " << scan.GetLevenbergLambda(); std::cout << ", original mss = " << mss << ", after update mss = " << new_mss << std::endl; std::cout.flush(); } // Set back parameters to old values // Set back parameters to old values for (unsigned int i=0; i& derivs, const NEWIMAGE::volume& mask, const NEWIMAGE::volume& mios, const NEWIMAGE::volume& pios, const NEWIMAGE::volume& jac, std::shared_ptr > susc, std::shared_ptr > bias, const NEWIMAGE::volume& pred, const NEWIMAGE::volume& dima, const NEWIMAGE::volume& sims, const NEWIMAGE::volume& pmask, const NEWMAT::Matrix& XtX, const NEWMAT::ColumnVector& Xty, const NEWMAT::ColumnVector& update) EddyTry { char fname[256], bname[256]; if (scan.IsSliceToVol()) strcpy(bname,"EDDY_DEBUG_S2V"); else strcpy(bname,"EDDY_DEBUG"); if (level>0) { sprintf(fname,"%s_masked_dima_%02d_%04d",bname,iter,scindx); NEWIMAGE::write_volume(dima*mask,fname); sprintf(fname,"%s_reverse_dima_%02d_%04d",bname,iter,scindx); NEWIMAGE::write_volume(pred-sims,fname); } if (level>1) { sprintf(fname,"%s_mask_%02d_%04d",bname,iter,scindx); NEWIMAGE::write_volume(mask,fname); sprintf(fname,"%s_pios_%02d_%04d",bname,iter,scindx); NEWIMAGE::write_volume(pios,fname); sprintf(fname,"%s_pred_%02d_%04d",bname,iter,scindx); NEWIMAGE::write_volume(pred,fname); sprintf(fname,"%s_dima_%02d_%04d",bname,iter,scindx); NEWIMAGE::write_volume(dima,fname); sprintf(fname,"%s_jac_%02d_%04d",bname,iter,scindx); NEWIMAGE::write_volume(jac,fname); sprintf(fname,"%s_orig_%02d_%04d",bname,iter,scindx); NEWIMAGE::write_volume(scan.GetIma(),fname); } if (level>2) { sprintf(fname,"%s_mios_%02d_%04d",bname,iter,scindx); NEWIMAGE::write_volume(mios,fname); sprintf(fname,"%s_pmask_%02d_%04d",bname,iter,scindx); NEWIMAGE::write_volume(pmask,fname); sprintf(fname,"%s_derivs_%02d_%04d",bname,iter,scindx); NEWIMAGE::write_volume(derivs,fname); sprintf(fname,"%s_XtX_%02d_%04d.txt",bname,iter,scindx); MISCMATHS::write_ascii_matrix(fname,XtX); sprintf(fname,"%s_Xty_%02d_%04d.txt",bname,iter,scindx); MISCMATHS::write_ascii_matrix(fname,Xty); sprintf(fname,"%s_update_%02d_%04d.txt",bname,iter,scindx); MISCMATHS::write_ascii_matrix(fname,update); } return; } EddyCatch std::vector EddyUtils::long_ec_update(// Input const std::vector >& pred, // Predictions in model space const NEWIMAGE::volume& pmask, // "Data valid" mask in model space float fwhm, // FWHM for Gaussian smoothing unsigned int nthr, // Number of threads to run it on bool very_verbose, // Write detailed output to screen? // These input parameters are for debugging only unsigned int iter, // Iteration unsigned int level, // Determines how much gets written const std::vector& debug_index, // Indicies of scans to write debug info for // Input/output EDDY::ECScanManager& sm) EddyTry // Scans we want to register to predictions { std::vector XtX(sm.NScans()); std::vector Xty(sm.NScans()); std::vector mss(sm.NScans(),0.0); std::vector updates; arma::mat joint_XtX; arma::colvec joint_Xty; // Divide the work of calculating the matrices for the updates if (very_verbose) cout << "Calculating updates for long time-constant EC" << endl; std::vector nvols = EddyUtils::ScansPerThread(sm.NScans(),nthr); std::vector threads(nthr-1); for (unsigned int i=0; i= total_old_mss) && sm.LongTimeConstantECModel().GetLevenbergLambda(0) < 1e6) { if (very_verbose) { // Write info about failed update if output is very verbose std::cout << "EddyUtils::long_ec_update: updates rejected for all volumes" << std::endl; std::cout << "EddyUtils::long_ec_update: Levenberg-lambda = " << sm.LongTimeConstantECModel().GetLevenbergLambda(0); std::cout << ", original mss = " << total_old_mss << ", after update mss = " << total_new_mss << std::endl; std::cout.flush(); } sm.LongTimeConstantECModel().UpdateParameters(sm,EddyUtils::Negate(updates)); // Set back to old values sm.LongTimeConstantECModel().SetLevenbergLambda(0,10.0*sm.LongTimeConstantECModel().GetLevenbergLambda(0)); regm = sm.LongTimeConstantECModel().GetLevenbergLambda(0) * arma::eye(nd,nd); // Tikhonov/Levenberg regularisation updates[0] = -arma::solve(joint_XtX+regm,joint_Xty,arma::solve_opts::likely_sympd); sm.LongTimeConstantECModel().UpdateParameters(sm,updates); total_new_mss = EddyUtils::calculate_long_ec_total_new_mss(pred,pmask,fwhm,sm,nthr,iter,level,debug_index); } if (std::isnan(total_new_mss) || total_new_mss >= total_old_mss) { // Set back if it never improved if (very_verbose) { // Write info about failed update if output is very verbose std::cout << "EddyUtils::long_ec_update: updates rejected for all volumes" << std::endl; std::cout << "EddyUtils::long_ec_update: Levenberg-lambda = " << sm.LongTimeConstantECModel().GetLevenbergLambda(0); std::cout << ", original mss = " << total_old_mss << ", after update mss = " << total_new_mss << std::endl; std::cout.flush(); } sm.LongTimeConstantECModel().UpdateParameters(sm,EddyUtils::Negate(updates)); // Set back to old values } } } else throw EddyException("EddyUtils::long_ec_update: Invalid long EC model"); return(mss); } EddyCatch void EddyUtils::long_ec_levenberg_updates(// Input unsigned int first_vol, // Index of first vol to work on unsigned int last_vol, // One past index of last vol to work on const std::vector >& pred, // Predictions in model space const NEWIMAGE::volume& pmask, // "Data valid" mask in model space float fwhm, const std::vector& XtX, // Hessian for esimation of updates const std::vector& Xty, // Gradient for estimation of updates const std::vector& mss, // mss for current long EC parameters bool very_verbose, // Should we write info to screen? // Input used for debugging output only unsigned int iter, // Iteration unsigned int level, // Determines how much gets written const std::vector& dbg_indx, // Indicies of scans to write debug info for // Input/output EDDY::ECScanManager& sm, // Scan manager // Output std::vector& updates) EddyTry// Applied updates { static std::mutex cout_mutex; // Magic static unsigned int nd = sm.LongTimeConstantECModel().NDerivs(); for (unsigned int s=first_vol; s= mss[s]) && sm.LongTimeConstantECModel().GetLevenbergLambda(s) < 1e6) { if (very_verbose) { // Write info about failed update if output is very verbose const std::lock_guard lock(cout_mutex); // Grab ownership of terminal output std::cout << "EddyUtils::long_ec_levenberg_updates: updates rejected for volume #" << s << std::endl; std::cout << "EddyUtils::long_ec_levenberg_updates: Levenberg-lambda = " << sm.LongTimeConstantECModel().GetLevenbergLambda(s); std::cout << ", original mss = " << mss[s] << ", after update mss = " << new_mss << std::endl; std::cout.flush(); } sm.LongTimeConstantECModel().UpdateParameters(sm.Scan(s),s,-updates[s]); // Set back to old values sm.LongTimeConstantECModel().SetLevenbergLambda(s,10.0*sm.LongTimeConstantECModel().GetLevenbergLambda(s)); // Increase lambda regm = sm.LongTimeConstantECModel().GetLevenbergLambda(s) * arma::eye(nd,nd); // Tikhonov/Levenberg regularisation updates[s] = -arma::solve(XtX[s]+regm,Xty[s],arma::solve_opts::likely_sympd); sm.LongTimeConstantECModel().UpdateParameters(sm.Scan(s),s,updates[s]); // Try the new updates new_mss = EddyUtils::calculate_new_mss(pred[s],sm.GetSuscHzOffResField(s),pmask,sm.Scan(s),fwhm,ParametersType::LongEC,s,iter,level,dbg_indx); } if (std::isnan(new_mss) || new_mss >= mss[s]) { // It means we failed if (very_verbose) { // Write info about failed update if output is very verbose const std::lock_guard lock(cout_mutex); // Grab ownership of terminal output std::cout << "EddyUtils::long_ec_levenberg_updates: updates rejected for volume #" << s << std::endl; std::cout << "EddyUtils::long_ec_levenberg_updates: Levenberg-lambda = " << sm.LongTimeConstantECModel().GetLevenbergLambda(s); std::cout << ", original mss = " << mss[s] << ", after update mss = " << new_mss << std::endl; std::cout.flush(); } sm.LongTimeConstantECModel().UpdateParameters(sm.Scan(s),s,-updates[s]); // Set back since it never improved updates[s].zeros(); // We did zero updates } } } return; } EddyCatch double EddyUtils::calculate_long_ec_total_new_mss(const std::vector >& pred, // Predictions in model space const NEWIMAGE::volume& pmask, // "Data valid" mask in model space float fwhm, // FWHM for Gaussian smoothing const EDDY::ECScanManager& sm, // Scan manager unsigned int nthr, // No. of threads // These input parameters are for debugging only unsigned int iter, // Iteration unsigned int level, // Determines how much gets written const std::vector& debug_index) EddyTry // Indicies of scans to write debug info for { std::vector mss(sm.NScans(),0.0); std::vector nvols = EddyUtils::ScansPerThread(sm.NScans(),nthr); std::vector threads(nthr-1); for (unsigned int i=0; i >& pred, // Predictions in model space const NEWIMAGE::volume& pmask, // "Data valid" mask in model space float fwhm, // FWHM for Gaussian smoothing const EDDY::ECScanManager& sm, // Scan manager // These input parameters are for debugging only unsigned int iter, // Iteration unsigned int level, // Determines how much gets written const std::vector& debug_indx, // Indicies of scans to write debug info for // Output std::vector mss) EddyTry // Vector of mean-sum-of-square-differences { for (unsigned int s=first_vol; s& pred, // Predictions in model space std::shared_ptr > susc, // Susceptibility-induced field const NEWIMAGE::volume& pmask, // "Data valid" mask in model space const EDDY::ECScan& scan, // Scan float fwhm, // FWHM for Gaussian smoothing EDDY::ParametersType whichp, // What parameters have been updated // These input parameters are for debugging only unsigned int scindx, // Scan index unsigned int iter, // Iteration unsigned int level, // Determines how much gets written const std::vector& debug_indx) EddyTry // Indicies of scans to write debug info for { NEWIMAGE::volume mask = pred; mask.setextrapolationmethod(NEWIMAGE::zeropad); mask = 0.0; NEWIMAGE::volume jac = pred; jac = 1.0; NEWIMAGE::volume pios = EddyUtils::transform_model_to_scan_space(pred,scan,susc,true,mask,&jac,NULL); // Transform binary mask into observation space NEWIMAGE::volume skrutt = pred; skrutt = 0.0; NEWIMAGE::volume mios = EddyUtils::transform_model_to_scan_space(pmask,scan,susc,false,skrutt,NULL,NULL); mios.binarise(0.99); // Value above (arbitrary) 0.99 implies valid voxels mask *= mios; // Volume and prediction mask falls within FOV // Calculate difference image between observed and predicted NEWIMAGE::volume dima = pios-scan.GetIma(); if (fwhm) { mask.setextrapolationmethod(NEWIMAGE::zeropad); dima = EddyUtils::Smooth(dima,fwhm,mask); } // Calculate mean-sum-of-squares from difference image double masksum = mask.sum(); double mss = (dima*mask).sumsquares() / masksum; if (level) { // If this is a debug run EddyUtils::write_debug_info_updated_mss(scan,scindx,whichp,debug_indx,iter,level,pred,mask,dima,pios,susc); } return(mss); } EddyCatch void EddyUtils::write_debug_info_updated_mss(const EDDY::ECScan& scan, unsigned int scindx, EDDY::ParametersType whichp, const std::vector& debug_indx, unsigned int iter, unsigned int level, const NEWIMAGE::volume& pred, const NEWIMAGE::volume& mask, const NEWIMAGE::volume& dima, const NEWIMAGE::volume& pios, std::shared_ptr > susc) EddyTry { if (std::any_of(debug_indx.begin(),debug_indx.end(),[scindx](unsigned int dindx){ return(dindx==scindx); })) { // If scindx among debug_indx if (level > 1) { char fname[256]; char const *basename = nullptr; if (whichp == EDDY::ParametersType::LongEC) basename = "EDDY_DEBUG_LONG_EC"; else { if (scan.IsSliceToVol()) basename = "EDDY_DEBUG_S2V"; else basename = "EDDY_DEBUG"; } sprintf(fname,"%s_new_masked_dima_%02d_%04d",basename,iter,scindx); NEWIMAGE::write_volume(dima*mask,fname); NEWIMAGE::volume sims = scan.GetUnwarpedIma(susc); sprintf(fname,"%s_new_reverse_dima_%02d_%04d",basename,iter,scindx); NEWIMAGE::write_volume(pred-sims,fname); sprintf(fname,"%s_new_mask_%02d_%04d",basename,iter,scindx); NEWIMAGE::write_volume(mask,fname); sprintf(fname,"%s_new_pios_%02d_%04d",basename,iter,scindx); NEWIMAGE::write_volume(pios,fname); sprintf(fname,"%s_new_dima_%02d_%04d",basename,iter,scindx); NEWIMAGE::write_volume(dima,fname); } } } EddyCatch void EddyUtils::long_ec_calculate_matrices(// Input unsigned int first_vol, // Index of first vol to work on unsigned int last_vol, // One past index of last vol to work on const std::vector >& pred, // Predictions in model space const NEWIMAGE::volume& pmask, // "Data valid" mask in model space float fwhm, // FWHM for Gaussian smoothing const EDDY::ECScanManager& sm, // Scan manager // These input parameters are for debugging only unsigned int iter, // Iteration unsigned int level, // Determines how much gets written const std::vector& debug_indx, // Indicies of scans to write debug info for // Output std::vector& mss, // Vector of sum-of-mean-square of difference std::vector& XtX, // Vector of XtX matrices std::vector& Xty) EddyTry // Vector of Xty vectors { for (unsigned int s=first_vol; s mask = pred[s]; mask.setextrapolationmethod(NEWIMAGE::zeropad); mask = 0.0; NEWIMAGE::volume jac = pred[s]; jac = 1.0; NEWIMAGE::volume pios = EddyUtils::transform_model_to_scan_space(pred[s],sm.Scan(s),sm.GetSuscHzOffResField(s),true,mask,&jac,NULL); // Transform binary mask into observation space NEWIMAGE::volume skrutt = pred[s]; skrutt = 0.0; NEWIMAGE::volume mios = EddyUtils::transform_model_to_scan_space(pmask,sm.Scan(s),sm.GetSuscHzOffResField(s),false,skrutt,NULL,NULL); mios.binarise(0.99); // Value above (arbitrary) 0.99 implies valid voxels mask *= mios; // Volume and prediction mask falls within FOV // Get partial derivatives w.r.t. to weights for previous and current EC field. NEWIMAGE::volume4D derivs = EddyUtils::get_long_ec_partial_derivatives_in_scan_space(pred[s],sm.Scan(s),s,sm.LongTimeConstantECModel(),sm.GetSuscHzOffResField(s)); // Smooth derivatives if needed if (fwhm) { mask.setextrapolationmethod(NEWIMAGE::zeropad); derivs = EddyUtils::Smooth(derivs,fwhm,mask); } // Calculate XtX where X is a matrix whose columns are the partial derivatives XtX[s] = EddyUtils::make_XtX(derivs,mask); // Calculate difference image between observed and predicted NEWIMAGE::volume dima = pios-sm.Scan(s).GetIma(); if (fwhm) { mask.setextrapolationmethod(NEWIMAGE::zeropad); dima = EddyUtils::Smooth(dima,fwhm,mask); } // Calculate Xty where y is the difference between observed and predicted. X as above. Xty[s] = EddyUtils::make_Xty(derivs,dima,mask); // Calculate update, only during testing period. Assuming individual weights arma::colvec update = -arma::solve(XtX[s],Xty[s],arma::solve_opts::likely_sympd); cout << "Previous = "; for (int ii=0; ii<22; ii++) cout << update(ii) << ", " << endl; cout << "Current = "; for (int ii=22; ii<44; ii++) cout << update(ii) << ", " << endl << std::flush; // Calculate mean-sum-of-squares from difference image double masksum = mask.sum(); mss[s] = (dima*mask).sumsquares() / masksum; // Write debug info if requested if (level) { NEWIMAGE::volume sims; // Used only for debug information sims = sm.Scan(s).GetUnwarpedIma(sm.GetSuscHzOffResField(s)); EddyUtils::write_debug_info_long_ec_calculate_matrices(sm.Scan(s),s,debug_indx,iter,level,derivs,mask,mios,pios, pred[s],dima,sims,pmask,XtX[s],Xty[s]); } } return; } EddyCatch void EddyUtils::write_debug_info_long_ec_calculate_matrices(const EDDY::ECScan& scan, unsigned int scindx, const std::vector& debug_indx, unsigned int iter, unsigned int level, const NEWIMAGE::volume4D& derivs, const NEWIMAGE::volume& mask, const NEWIMAGE::volume& mios, const NEWIMAGE::volume& pios, const NEWIMAGE::volume& pred, const NEWIMAGE::volume& dima, const NEWIMAGE::volume& sims, const NEWIMAGE::volume& pmask, const arma::mat& XtX, const arma::colvec& Xty) EddyTry { if (std::any_of(debug_indx.begin(),debug_indx.end(),[scindx](unsigned int dindx){ return(dindx==scindx); })) { // If scindx among debug_indx char fname[256]; if (level > 0) { sprintf(fname,"EDDY_DEBUG_LONG_EC_masked_dima_%02d_%04d",iter,scindx); NEWIMAGE::write_volume(dima*mask,fname); sprintf(fname,"EDDY_DEBUG_LONG_EC_reverse_dima_%02d_%04d",iter,scindx); NEWIMAGE::write_volume(pred-sims,fname); } if (level > 1) { sprintf(fname,"EDDY_DEBUG_LONG_EC_mask_%02d_%04d",iter,scindx); NEWIMAGE::write_volume(mask,fname); sprintf(fname,"EDDY_DEBUG_LONG_EC_pios_%02d_%04d",iter,scindx); NEWIMAGE::write_volume(pios,fname); sprintf(fname,"EDDY_DEBUG_LONG_EC_pred_%02d_%04d",iter,scindx); NEWIMAGE::write_volume(pred,fname); sprintf(fname,"EDDY_DEBUG_LONG_EC_dima_%02d_%04d",iter,scindx); NEWIMAGE::write_volume(dima,fname); sprintf(fname,"EDDY_DEBUG_LONG_EC_orig_%02d_%04d",iter,scindx); NEWIMAGE::write_volume(scan.GetIma(),fname); } if (level > 2) { sprintf(fname,"EDDY_DEBUG_LONG_EC_mios_%02d_%04d",iter,scindx); NEWIMAGE::write_volume(mios,fname); sprintf(fname,"EDDY_DEBUG_LONG_EC_pmask_%02d_%04d",iter,scindx); NEWIMAGE::write_volume(pmask,fname); sprintf(fname,"EDDY_DEBUG_LONG_EC_derivs_%02d_%04d",iter,scindx); NEWIMAGE::write_volume(derivs,fname); sprintf(fname,"EDDY_DEBUG_LONG_EC_XtX_%02d_%04d",iter,scindx); MISCMATHS::write_ascii_matrix(fname,XtX); sprintf(fname,"EDDY_DEBUG_LONG_EC_Xty_%02d_%04d",iter,scindx); MISCMATHS::write_ascii_matrix(fname,Xty); } } } EddyCatch /* double EddyUtils::param_update(// Input const NEWIMAGE::volume& pred, // Prediction in model space boost::shared_ptr > susc, // Susceptibility induced off-resonance field const NEWIMAGE::volume& pmask, // Parameters whichp, // What parameters do we want to update bool cbs, // Check (success of parameters) Before Set // Input/output EDDY::ECScan& scan, // Scan we want to register to pred // Output NEWMAT::ColumnVector *rupdate) // Vector of updates, optional output { // Transform prediction into observation space NEWIMAGE::volume pios = EddyUtils::TransformModelToScanSpace(pred,scan,susc); // Transform binary mask into observation space NEWIMAGE::volume mask = pred; mask = 0.0; NEWIMAGE::volume bios = EddyUtils::transform_model_to_scan_space(pmask,scan,susc,false,mask,NULL,NULL); bios.binarise(0.99); // Value above (arbitrary) 0.99 implies valid voxels mask *= bios; // Volume and prediction mask falls within FOV // Get partial derivatives w.r.t. requested category of parameters in prediction space NEWIMAGE::volume4D derivs = EddyUtils::get_partial_derivatives_in_scan_space(pred,scan,susc,whichp); // Calculate XtX where X is a matrix whos columns are the partial derivatives NEWMAT::Matrix XtX = EddyUtils::make_XtX(derivs,mask); // Calculate difference image between observed and predicted NEWIMAGE::volume dima = pios-scan.GetIma(); // Calculate Xty where y is the difference between observed and predicted. X as above. NEWMAT::ColumnVector Xty = EddyUtils::make_Xty(derivs,dima,mask); // Calculate mean sum of squares from difference image double mss = (dima*mask).sumsquares() / mask.sum(); // Calculate update to parameters NEWMAT::ColumnVector update = -XtX.i()*Xty; // Update parameters for (unsigned int i=0; i mss) { // Oh dear for (unsigned int i=0; i EddyUtils::transform_model_to_scan_space(// Input const NEWIMAGE::volume& pred, const EDDY::ECScan& scan, std::shared_ptr > susc, bool jacmod, // Output NEWIMAGE::volume& omask, NEWIMAGE::volume *jac, NEWIMAGE::volume4D *grad) EddyTry { // Get total field from scan if (jacmod && !jac) throw EddyException("EddyUtils::transform_model_to_scan_space: jacmod can only be used with valid jac"); NEWIMAGE::volume4D dfield; if (jacmod || jac) dfield = scan.FieldForModelToScanTransform(susc,omask,*jac); else dfield = scan.FieldForModelToScanTransform(susc,omask); NEWMAT::Matrix eye(4,4); eye=0; eye(1,1)=1.0; eye(2,2)=1.0; eye(3,3)=1.0; eye(4,4)=1.0; NEWIMAGE::volume ovol = pred; ovol = 0.0; NEWIMAGE::volume mask3(pred.xsize(),pred.ysize(),pred.zsize()); NEWIMAGE::copybasicproperties(pred,mask3); mask3 = 0; std::vector ddir(3); ddir[0] = 0; ddir[1] = 1; ddir[2] = 2; if (scan.IsSliceToVol()) { if (grad) { grad->reinitialize(pred.xsize(),pred.ysize(),pred.zsize(),3); NEWIMAGE::copybasicproperties(pred,*grad); } for (unsigned int tp=0; tp slices = scan.GetMBG().SlicesAtTimePoint(tp); if (grad) NEWIMAGE::raw_general_transform(pred,eye,dfield,ddir,ddir,slices,&eye,&R,ovol,*grad,&mask3); else NEWIMAGE::apply_warp(pred,eye,dfield,eye,R,slices,ovol,mask3); } } else { std::vector all_slices; // Get RB matrix NEWMAT::Matrix R = scan.ForwardMovementMatrix(); // Transform prediction using RB, inverted Tot map and Jacobian if (grad) { grad->reinitialize(pred.xsize(),pred.ysize(),pred.zsize(),3); NEWIMAGE::copybasicproperties(pred,*grad); NEWIMAGE::raw_general_transform(pred,eye,dfield,ddir,ddir,all_slices,&eye,&R,ovol,*grad,&mask3); } else NEWIMAGE::apply_warp(pred,eye,dfield,eye,R,ovol,mask3); } omask *= EddyUtils::ConvertMaskToFloat(mask3); // Combine all masks EddyUtils::SetTrilinearInterp(omask); if (jacmod) ovol *= *jac; // Jacobian modulation if it was asked for return(ovol); } EddyCatch // This is a temporary version to allow for writing of debug information /* NEWIMAGE::volume EddyUtils::transform_model_to_scan_space(// Input const NEWIMAGE::volume& pred, const EDDY::ECScan& scan, std::shared_ptr > susc, bool jacmod, // Output NEWIMAGE::volume& omask, NEWIMAGE::volume *jac, NEWIMAGE::volume4D *grad, // Tmp unsigned int scindx, unsigned int iter, unsigned int level) EddyTry { // Get total field from scan if (jacmod && !jac) throw EddyException("EddyUtils::transform_model_to_scan_space: jacmod can only be used with valid jac"); NEWIMAGE::volume4D dfield; if (jacmod || jac) dfield = scan.FieldForModelToScanTransform(susc,omask,*jac); else dfield = scan.FieldForModelToScanTransform(susc,omask); NEWMAT::Matrix eye(4,4); eye=0; eye(1,1)=1.0; eye(2,2)=1.0; eye(3,3)=1.0; eye(4,4)=1.0; NEWIMAGE::volume ovol = pred; ovol = 0.0; NEWIMAGE::volume mask3(pred.xsize(),pred.ysize(),pred.zsize()); NEWIMAGE::copybasicproperties(pred,mask3); mask3 = 0; std::vector ddir(3); ddir[0] = 0; ddir[1] = 1; ddir[2] = 2; if (scan.IsSliceToVol()) { if (grad) { grad->reinitialize(pred.xsize(),pred.ysize(),pred.zsize(),3); NEWIMAGE::copybasicproperties(pred,*grad); } for (unsigned int tp=0; tp slices = scan.GetMBG().SlicesAtTimePoint(tp); if (grad) NEWIMAGE::raw_general_transform(pred,eye,dfield,ddir,ddir,slices,&eye,&R,ovol,*grad,&mask3); else NEWIMAGE::apply_warp(pred,eye,dfield,eye,R,slices,ovol,mask3); } } else { std::vector all_slices; // Get RB matrix NEWMAT::Matrix R = scan.ForwardMovementMatrix(); // Transform prediction using RB, inverted Tot map and Jacobian if (grad) { grad->reinitialize(pred.xsize(),pred.ysize(),pred.zsize(),3); NEWIMAGE::copybasicproperties(pred,*grad); NEWIMAGE::raw_general_transform(pred,eye,dfield,ddir,ddir,all_slices,&eye,&R,ovol,*grad,&mask3); } else NEWIMAGE::apply_warp(pred,eye,dfield,eye,R,ovol,mask3); } // mask3 is buggered already here char bfname[256]; char fname[256]; if (level) { if (scan.IsSliceToVol()) strcpy(bfname,"EDDY_DEBUG_SPECIAL_S2V"); else strcpy(bfname,"EDDY_DEBUG_SPECIAL"); sprintf(fname,"%s_omask_before_%02d_%04d",bfname,iter,scindx); NEWIMAGE::write_volume(omask,fname); sprintf(fname,"%s_mask3_before_%02d_%04d",bfname,iter,scindx); NEWIMAGE::write_volume(mask3,fname); } omask *= EddyUtils::ConvertMaskToFloat(mask3); // Combine all masks if (level) { sprintf(fname,"%s_omask_after_1_%02d_%04d",bfname,iter,scindx); NEWIMAGE::write_volume(omask,fname); } EddyUtils::SetTrilinearInterp(omask); if (level) { sprintf(fname,"%s_omask_after_2_%02d_%04d",bfname,iter,scindx); NEWIMAGE::write_volume(omask,fname); } if (jacmod) ovol *= *jac; // Jacobian modulation if it was asked for return(ovol); } EddyCatch */ // Has been modified for slice-to-vol EDDY::ImageCoordinates EddyUtils::transform_coordinates_from_model_to_scan_space(// Input const NEWIMAGE::volume& pred, const EDDY::ECScan& scan, std::shared_ptr > susc, // Output NEWIMAGE::volume *omask, NEWIMAGE::volume *jac) EddyTry { // Get total field from scan NEWIMAGE::volume4D dfield; if (omask && jac) dfield = scan.FieldForModelToScanTransform(susc,*omask,*jac); else if (omask) dfield = scan.FieldForModelToScanTransform(susc,*omask); else if (jac) dfield = scan.FieldForModelToScanTransformWithJac(susc,*jac); else dfield = scan.FieldForModelToScanTransform(susc); ImageCoordinates coord(pred.xsize(),pred.ysize(),pred.zsize()); if (scan.IsSliceToVol()) { for (unsigned int tp=0; tp slices = scan.GetMBG().SlicesAtTimePoint(tp); // Transform coordinates using RB and inverted Tot map if (omask) { NEWIMAGE::volume mask2(pred.xsize(),pred.ysize(),pred.zsize()); NEWIMAGE::copybasicproperties(pred,mask2); mask2 = 0; EddyUtils::transform_coordinates(pred,dfield,R,slices,coord,&mask2); *omask *= mask2; EddyUtils::SetTrilinearInterp(*omask); } else EddyUtils::transform_coordinates(pred,dfield,R,slices,coord,NULL); } } else { // Get RB matrix NEWMAT::Matrix R = scan.ForwardMovementMatrix(); std::vector all_slices; // Transform coordinates using RB and inverted Tot map if (omask) { NEWIMAGE::volume mask2(pred.xsize(),pred.ysize(),pred.zsize()); NEWIMAGE::copybasicproperties(pred,mask2); mask2 = 0; EddyUtils::transform_coordinates(pred,dfield,R,all_slices,coord,&mask2); *omask *= mask2; EddyUtils::SetTrilinearInterp(*omask); } else EddyUtils::transform_coordinates(pred,dfield,R,all_slices,coord,NULL); } return(coord); } EddyCatch // Has been modified for slice-to-vol NEWIMAGE::volume4D EddyUtils::get_partial_derivatives_in_scan_space(// Input const NEWIMAGE::volume& pred, // Prediction in model space const EDDY::ECScan& scan, // Scan space std::shared_ptr > susc, // Susceptibility off-resonance field EDDY::ParametersType whichp) EddyTry { NEWIMAGE::volume basejac; NEWIMAGE::volume4D grad; NEWIMAGE::volume skrutt(pred.xsize(),pred.ysize(),pred.zsize()); NEWIMAGE::volume base = transform_model_to_scan_space(pred,scan,susc,true,skrutt,&basejac,&grad); ImageCoordinates basecoord = transform_coordinates_from_model_to_scan_space(pred,scan,susc,NULL,NULL); NEWIMAGE::volume4D derivs(base.xsize(),base.ysize(),base.zsize(),scan.NDerivs(whichp)); NEWIMAGE::copybasicproperties(scan.GetIma(),derivs); NEWIMAGE::volume jac = pred; ECScan sc = scan; // We are relying on the order of derivatives being movement followed by EC. if (whichp == EDDY::ParametersType::All || whichp == EDDY::ParametersType::Movement) { // If we are asked for movement derivatives // First we calculate the movement derivatives using modulation. for (unsigned int i=0; i smod = di.GetSliceModulator(j).GetMod(); for (int sl=0; sl EddyUtils::get_partial_derivatives_in_scan_space(// Input const NEWIMAGE::volume& pred, // Prediction in model space const EDDY::ECScan& scan, // Scan space std::shared_ptr > susc, // Susceptibility off-resonance field EDDY::Parameters whichp) { NEWIMAGE::volume basejac; NEWIMAGE::volume4D grad; NEWIMAGE::volume skrutt(pred.xsize(),pred.ysize(),pred.zsize()); NEWIMAGE::volume base = transform_model_to_scan_space(pred,scan,susc,true,skrutt,&basejac,&grad); ImageCoordinates basecoord = transform_coordinates_from_model_to_scan_space(pred,scan,susc,NULL,NULL); NEWIMAGE::volume4D derivs(base.xsize(),base.ysize(),base.zsize(),scan.NDerivs(whichp)); NEWIMAGE::volume jac = pred; ECScan sc = scan; for (unsigned int i=0; i EddyUtils::get_direct_partial_derivatives_in_scan_space(// Input const NEWIMAGE::volume& pred, // Prediction in model space const EDDY::ECScan& scan, // Scan space std::shared_ptr > susc, // Susceptibility off-resonance field EDDY::ParametersType whichp) EddyTry { NEWIMAGE::volume jac(pred.xsize(),pred.ysize(),pred.zsize()); NEWIMAGE::volume skrutt(pred.xsize(),pred.ysize(),pred.zsize()); NEWIMAGE::volume base = transform_model_to_scan_space(pred,scan,susc,true,skrutt,&jac,NULL); NEWIMAGE::volume4D derivs(base.xsize(),base.ysize(),base.zsize(),scan.NDerivs(whichp)); ECScan sc = scan; for (unsigned int i=0; i perturbed = transform_model_to_scan_space(pred,sc,susc,true,skrutt,&jac,NULL); derivs[i] = (perturbed-base) / sc.GetDerivScale(i,whichp); sc.SetDerivParam(i,p,whichp); } return(derivs); } EddyCatch NEWIMAGE::volume4D EddyUtils::get_long_ec_partial_derivatives_in_scan_space(// Input const NEWIMAGE::volume& pred, // Prediction in model space const EDDY::ECScan& scan, // Scan space unsigned int si, // Scan index const EDDY::LongECModel& lecm, // Model for long EC std::shared_ptr > susc) // Susceptibility off-resonance field { NEWIMAGE::volume basejac; // Jacobian map with input parameters NEWIMAGE::volume4D grad; // Gradient with input parameters NEWIMAGE::volume skrutt(pred.xsize(),pred.ysize(),pred.zsize()); // pred in scan-space with input parameters NEWIMAGE::volume base = transform_model_to_scan_space(pred,scan,susc,true,skrutt,&basejac,&grad); // Sampling points with input parameters ImageCoordinates basecoord = transform_coordinates_from_model_to_scan_space(pred,scan,susc,NULL,NULL); NEWIMAGE::volume jac = pred; // Jacobian map with the different perturbations ECScan sc = scan; // Copy to allow for changes to parameters // The derivatives we want to calculate NEWIMAGE::volume4D derivs(base.xsize(),base.ysize(),base.zsize(),lecm.NDerivs()); for (unsigned int i=0; i EddyUtils::get_TV_EC_weights_partial_derivatives_in_scan_space(// Input const NEWIMAGE::volume& pred, // Prediction in model space const EDDY::ECScan& scan, // Scan space std::shared_ptr > susc) // Susceptibility off-resonance field { const double perturbation = 1e-6; // Should probably be a property of ECScan/ScanManager instead. NEWIMAGE::volume basejac; // Jacobian map with input parameters NEWIMAGE::volume4D grad; // Gradient with input parameters NEWIMAGE::volume skrutt(pred.xsize(),pred.ysize(),pred.zsize()); // pred in scan-space with input parameters NEWIMAGE::volume base = transform_model_to_scan_space(pred,scan,susc,true,skrutt,&basejac,&grad); // Sampling points with input parameters ImageCoordinates basecoord = transform_coordinates_from_model_to_scan_space(pred,scan,susc,NULL,NULL); NEWIMAGE::volume jac = pred; // Jacobian map with the different perturbations ECScan sc = scan; // Copy to allow for changes to parameters // The derivatives we want to calculate NEWIMAGE::volume4D derivs(base.xsize(),base.ysize(),base.zsize(),2*sc.GetMBG().NGroups()); // Derivatives for weights of EC from previous volume std::vector weights = sc.GetWeightsForPrevious(); for (unsigned int i=0; i EddyUtils::get_TV_EC_time_constant_derivative_in_scan_space(// Input const NEWIMAGE::volume& pred, // Prediction in model space const EDDY::ECScan& scan, // Scan space std::shared_ptr > susc) // Susceptibility off-resonance field { const double perturbation = 1e-5; // Should probably be a property of ECScan/ScanManager instead. NEWIMAGE::volume basejac; // Jacobian map with input parameters NEWIMAGE::volume4D grad; // Gradient with input parameters NEWIMAGE::volume skrutt(pred.xsize(),pred.ysize(),pred.zsize()); // pred in scan-space with input parameters NEWIMAGE::volume base = transform_model_to_scan_space(pred,scan,susc,true,skrutt,&basejac,&grad); // Sampling points with input parameters ImageCoordinates basecoord = transform_coordinates_from_model_to_scan_space(pred,scan,susc,NULL,NULL); NEWIMAGE::volume jac = pred; // Jacobian map with the different perturbations // The derivative we want to calculate NEWIMAGE::volume4D deriv(base.xsize(),base.ysize(),base.zsize(),1); // Strictly speaking 3D ECScan sc = scan; // Copy to allow for changes to parameter // And finally calculate derivative sc.SetTimeConstantForLongEC(sc.GetTimeConstantForLongEC()+perturbation); ImageCoordinates diff = transform_coordinates_from_model_to_scan_space(pred,sc,susc,NULL,&jac) - basecoord; deriv[0] = (diff*grad) / perturbation; deriv[0] += base * (jac-basejac) / perturbation; return(deriv); } */ /* NEWIMAGE::volume EddyUtils::TransformScanToModelSpace(// Input const EDDY::ECScan& scan, boost::shared_ptr > susc, // Output NEWIMAGE::volume& omask) { // Get total field from scan NEWIMAGE::volume jac; NEWIMAGE::volume4D dfield = scan.FieldForScanToModelTransform(susc,omask,jac); // Transform prediction using inverse RB, dfield and Jacobian NEWMAT::Matrix iR = scan.InverseMovementMatrix(); NEWIMAGE::volume ovol = scan.GetIma(); ovol = 0.0; NEWIMAGE::volume mask2(ovol.xsize(),ovol.ysize(),ovol.zsize()); NEWIMAGE::copybasicproperties(scan.GetIma(),mask2); mask2 = 1; NEWIMAGE::general_transform(scan.GetIma(),iR,dfield,ovol,mask2); // Combine all masks omask *= EddyUtils::ConvertMaskToFloat(mask2); return(jac*ovol); } */ NEWIMAGE::volume EddyUtils::DirectTransformScanToModelSpace(// Input const EDDY::ECScan& scan, std::shared_ptr > susc, // Output NEWIMAGE::volume& omask) EddyTry { NEWIMAGE::volume ima = scan.GetIma(); NEWIMAGE::volume eb = scan.ECField(); NEWIMAGE::volume4D dfield = FieldUtils::Hz2VoxelDisplacements(eb,scan.GetAcqPara()); dfield = FieldUtils::Voxel2MMDisplacements(dfield); NEWMAT::Matrix eye(4,4); eye=0; eye(1,1)=1.0; eye(2,2)=1.0; eye(3,3)=1.0; eye(4,4)=1.0; NEWIMAGE::volume ovol = ima; ovol = 0.0; NEWIMAGE::volume mask(ima.xsize(),ima.ysize(),ima.zsize()); NEWIMAGE::apply_warp(ima,eye,dfield,eye,eye,ovol,mask); NEWMAT::Matrix iR = scan.InverseMovementMatrix(); NEWIMAGE::volume tmp = ovol; ovol = 0; NEWIMAGE::affine_transform(tmp,iR,ovol,mask); return(ovol); } EddyCatch // Right now the following function is clearly malformed. It is missing susc /* NEWIMAGE::volume EddyUtils::DirectTransformModelToScanSpace(// Input const NEWIMAGE::volume& ima, const EDDY::ECScan& scan, const EDDY::MultiBandGroups& mbg, boost::shared_ptr > susc, // Output NEWIMAGE::volume& omask) { if (scan.IsSliceToVol()) { NEWMAT::Matrix R = scan.ForwardMovementMatrix(); NEWIMAGE::volume tmp = ima; tmp = 0; NEWIMAGE::volume mask(tmp.xsize(),tmp.ysize(),tmp.zsize()); NEWIMAGE::affine_transform(ima,R,tmp,mask); } else { } NEWIMAGE::volume eb = scan.ECField(); NEWIMAGE::volume4D dfield = FieldUtils::Hz2VoxelDisplacements(eb,scan.GetAcqPara()); NEWIMAGE::volume4D idfield = FieldUtils::InvertDisplacementField(dfield,scan.GetAcqPara(),EddyUtils::ConvertMaskToFloat(mask),omask); idfield = FieldUtils::Voxel2MMDisplacements(idfield); NEWMAT::Matrix eye(4,4); eye=0; eye(1,1)=1.0; eye(2,2)=1.0; eye(3,3)=1.0; eye(4,4)=1.0; NEWIMAGE::volume ovol = tmp; ovol = 0.0; NEWIMAGE::apply_warp(tmp,eye,idfield,eye,eye,ovol,mask); return(ovol); } */ NEWMAT::Matrix EddyUtils::make_XtX(const NEWIMAGE::volume4D& vols, const NEWIMAGE::volume& mask) EddyTry { NEWMAT::Matrix XtX(vols.tsize(),vols.tsize()); XtX = 0.0; for (int r=1; r<=vols.tsize(); r++) { for (int c=r; c<=vols.tsize(); c++) { for (NEWIMAGE::volume::fast_const_iterator rit=vols.fbegin(r-1), ritend=vols.fend(r-1), cit=vols.fbegin(c-1), mit=mask.fbegin(); rit!=ritend; ++rit, ++cit, ++mit) { if (*mit) XtX(r,c) += (*rit)*(*cit); } } } for (int r=2; r<=vols.tsize(); r++) { for (int c=1; c& Xvols, const NEWIMAGE::volume& Yvol, const NEWIMAGE::volume& mask) EddyTry { NEWMAT::ColumnVector Xty(Xvols.tsize()); Xty = 0.0; for (int r=1; r<=Xvols.tsize(); r++) { for (NEWIMAGE::volume::fast_const_iterator Xit=Xvols.fbegin(r-1), Xend=Xvols.fend(r-1), Yit=Yvol.fbegin(), mit=mask.fbegin(); Xit!=Xend; ++Xit, ++Yit, ++mit) { if (*mit) Xty(r) += (*Xit)*(*Yit); } } return(Xty); } EddyCatch // Has been modified for slice-to-vol void EddyUtils::transform_coordinates(// Input const NEWIMAGE::volume& f, const NEWIMAGE::volume4D& d, const NEWMAT::Matrix& M, std::vector slices, // Input/Output ImageCoordinates& c, // Output NEWIMAGE::volume *omask) EddyTry { NEWMAT::Matrix iA = d[0].sampling_mat(); float A11=iA(1,1), A12=iA(1,2), A13=iA(1,3), A14=iA(1,4); float A21=iA(2,1), A22=iA(2,2), A23=iA(2,3), A24=iA(2,4); float A31=iA(3,1), A32=iA(3,2), A33=iA(3,3), A34=iA(3,4); // Create a matrix mapping from mm-coordinates in volume i // to voxel coordinates in volume f. If the matrix M is empty // this is simply a mm->voxel mapping for volume f NEWMAT::Matrix iM = f.sampling_mat().i() * M.i(); float M11=iM(1,1), M12=iM(1,2), M13=iM(1,3), M14=iM(1,4); float M21=iM(2,1), M22=iM(2,2), M23=iM(2,3), M24=iM(2,4); float M31=iM(3,1), M32=iM(3,2), M33=iM(3,3), M34=iM(3,4); // If no slices were specified, do all slices if (slices.size() == 0) { slices.resize(c.NZ()); for (unsigned int z=0; z c.NZ()) throw EddyException("EddyUtils::transform_coordinates: slices vector too long"); else { for (unsigned int z=0; z= c.NZ()) throw EddyException("EddyUtils::transform_coordinates: slices vector has invalid entry");} for (unsigned int k=0; k maxpetr) return(false); for (unsigned int i=morder; i<3*morder; i++) if (std::abs(update(i+1)) > maxtr) return(false); } else { for (unsigned int i=0; i maxtr) return(false); for (unsigned int i=morder; i<2*morder; i++) if (std::abs(update(i+1)) > maxpetr) return(false); for (unsigned int i=2*morder; i<3*morder; i++) if (std::abs(update(i+1)) > maxtr) return(false); } // Rotations for (unsigned int i=3*morder; i<6*morder; i++) { if (std::abs(update(i+1)) > maxrot) return(false); } // First order EC if (scan.Model() == ECModelType::Linear || scan.Model() == ECModelType::Quadratic || scan.Model() == ECModelType::Cubic) { for (unsigned int i=6*morder; i<6*morder+3; i++) { if (std::abs(update(i+1)) > maxfirst) return(false); } } // Second order EC if (scan.Model() == ECModelType::Quadratic || scan.Model() == ECModelType::Cubic) { for (unsigned int i=6*morder+3; i<6*morder+9; i++) { if (std::abs(update(i+1)) > maxsecond) return(false); } } // I currently lack data for 3rd order EC return(true); } EddyCatch //@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@ // // Class FieldUtils // // Helper Class used to perform various useful calculations // on displacement fields. // //@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@ NEWIMAGE::volume4D FieldUtils::Hz2VoxelDisplacements(const NEWIMAGE::volume& hzfield, const AcqPara& acqp) EddyTry { NEWIMAGE::volume4D dfield(hzfield.xsize(),hzfield.ysize(),hzfield.zsize(),3); NEWIMAGE::copybasicproperties(hzfield,dfield); for (int i=0; i<3; i++) dfield[i] = float((acqp.PhaseEncodeVector())(i+1) * acqp.ReadOutTime()) * hzfield; return(dfield); } EddyCatch NEWIMAGE::volume4D FieldUtils::Hz2MMDisplacements(const NEWIMAGE::volume& hzfield, const AcqPara& acqp) EddyTry { NEWIMAGE::volume4D dfield(hzfield.xsize(),hzfield.ysize(),hzfield.zsize(),3); NEWIMAGE::copybasicproperties(hzfield,dfield); dfield[0] = float(hzfield.xdim()*(acqp.PhaseEncodeVector())(1) * acqp.ReadOutTime()) * hzfield; dfield[1] = float(hzfield.ydim()*(acqp.PhaseEncodeVector())(2) * acqp.ReadOutTime()) * hzfield; dfield[2] = float(hzfield.zdim()*(acqp.PhaseEncodeVector())(3) * acqp.ReadOutTime()) * hzfield; return(dfield); } EddyCatch ///////////////////////////////////////////////////////////////////// // // Inverts a 1D displacementfield. The input field should be in units // of voxels and the output will be too. // ///////////////////////////////////////////////////////////////////// NEWIMAGE::volume FieldUtils::Invert1DDisplacementField(// Input const NEWIMAGE::volume& dfield, const AcqPara& acqp, const NEWIMAGE::volume& inmask, // Output NEWIMAGE::volume& omask) EddyTry { NEWIMAGE::volume fc = dfield; // fc : field copy NEWIMAGE::volume imc = inmask; // imc: inmask copy // Make it so that we invert in first (x) direction unsigned int d=0; for (; d<3; d++) if ((acqp.PhaseEncodeVector())(d+1)) break; if (d==1) { fc.swapdimensions(2,1,3); imc.swapdimensions(2,1,3); omask.swapdimensions(2,1,3); } else if (d==2) { fc.swapdimensions(3,2,1); imc.swapdimensions(3,2,1); omask.swapdimensions(3,2,1); } NEWIMAGE::volume idf = fc; // idf : inverse displacement field // Do the inversion for (int k=0; k0 && ii0; ii--) idf(ii-1,j,k) = idf(ii,j,k); // Process NaN's at end of column for (ii=idf.xsize()-1; ii>0 && idf(ii,j,k)==FLT_MAX; ii--) ; // On purpose for (; ii FieldUtils::Invert3DDisplacementField(// Input const NEWIMAGE::volume4D& dfield, const AcqPara& acqp, const NEWIMAGE::volume& inmask, // Output NEWIMAGE::volume& omask) EddyTry { NEWIMAGE::volume4D idfield = dfield; idfield = 0.0; unsigned int cnt=0; for (unsigned int i=0; i<3; i++) if ((acqp.PhaseEncodeVector())(i+1)) cnt++; if (cnt != 1) throw EddyException("FieldUtils::InvertDisplacementField: Phase encode vector must have exactly one non-zero component"); unsigned int i=0; for (; i<3; i++) if ((acqp.PhaseEncodeVector())(i+1)) break; idfield[i] = Invert1DDisplacementField(dfield[i],acqp,inmask,omask); return(idfield); } EddyCatch ///////////////////////////////////////////////////////////////////// // // Calculates the Jacobian determinant of a 3D displacement field. // The field must be in units of voxels and in the present // implementation it must also be inherently 1D. // ///////////////////////////////////////////////////////////////////// NEWIMAGE::volume FieldUtils::GetJacobian(const NEWIMAGE::volume4D& dfield, const AcqPara& acqp) EddyTry { unsigned int cnt=0; for (unsigned int i=0; i<3; i++) if ((acqp.PhaseEncodeVector())(i+1)) cnt++; if (cnt != 1) throw EddyException("FieldUtils::GetJacobian: Phase encode vector must have exactly one non-zero component"); unsigned int i=0; for (; i<3; i++) if ((acqp.PhaseEncodeVector())(i+1)) break; NEWIMAGE::volume jacfield = GetJacobianFrom1DField(dfield[i],i); return(jacfield); } EddyCatch ///////////////////////////////////////////////////////////////////// // // Calculates the Jacobian determinant of a 1D displacement field. // The field must be in units of voxels. // ///////////////////////////////////////////////////////////////////// NEWIMAGE::volume FieldUtils::GetJacobianFrom1DField(const NEWIMAGE::volume& dfield, unsigned int dir) EddyTry { // Calculate spline coefficients for displacement field std::vector dim(3,0); dim[0] = dfield.xsize(); dim[1] = dfield.ysize(); dim[2] = dfield.zsize(); std::vector ep(3,SPLINTERPOLATOR::Mirror); SPLINTERPOLATOR::Splinterpolator spc(dfield.fbegin(),dim,ep,3,false); // Get Jacobian at voxel centres NEWIMAGE::volume jacf = dfield; for (int k=0; k icsl; if (_sm.MultiBand().MBFactor() == 1) icsl = _sm.IntraCerebralSlices(500); // N.B. Hardcoded. Might want to move up _tr.ReSize(3,_sm.NScans(ScanType::Any)); _rot.ReSize(3,_sm.NScans(ScanType::Any)); for (unsigned int i=0; i<_sm.NScans(ScanType::Any); i++) { for (unsigned int j=0; j<3; j++) { _tr(j+1,i+1) = _sm.Scan(i,ScanType::Any).GetMovementStd(j,icsl); _rot(j+1,i+1) = 180.0 * _sm.Scan(i,ScanType::Any).GetMovementStd(3+j,icsl) / 3.141592653589793; } } } EddyCatch /****************************************************************//** * * \brief Returns a vector of indices to volumes with little movement * * * ********************************************************************/ std::vector s2vQuant::FindStillVolumes(ScanType st, const std::vector& mbsp) const EddyTry { std::vector rval; for (unsigned int i=0; i<_sm.NScans(st); i++) { unsigned int j = i; if (st==ScanType::b0) j = _sm.Getb02GlobalIndexMapping(i); else if (st==ScanType::DWI) j = _sm.GetDwi2GlobalIndexMapping(i); bool is_still = true; for (unsigned int pi=0; pi _trth) is_still = false; } else if (mbsp[pi] > 2 && mbsp[pi] < 6) { NEWMAT::ColumnVector tmp = _rot.Column(j+1); if (tmp(mbsp[pi]+1-3) > _rotth) is_still = false; } else throw EddyException("s2vQuant::FindStillVolumes: mbsp out of range"); } if (is_still) rval.push_back(i); } return(rval); } EddyCatch