int w2xconv_convert_rgb(struct W2XConv *conv, unsigned char *dst, size_t dst_step_byte, /* rgb24 (src_w*ratio, src_h*ratio) */ unsigned char *src, size_t src_step_byte, /* rgb24 (src_w, src_h) */ int src_w, int src_h, int denoise_level, /* 0:none, 1:L1 denoise, other:L2 denoise */ double scale, int block_size) { int dst_h = src_h * scale; int dst_w = src_w * scale; cv::Mat srci(src_h, src_w, CV_8UC3, src, src_step_byte); cv::Mat dsti(dst_h, dst_w, CV_8UC3, dst, dst_step_byte); cv::Mat image; bool is_rgb = (conv->impl->scale2_models[0]->getNInputPlanes() == 3); if (is_rgb) { srci.copyTo(image); convert_mat(conv, image, denoise_level, scale, dst_w, dst_h, block_size, w2xc::IMAGE_RGB); image.copyTo(dsti); } else { srci.convertTo(image, CV_32F, 1.0 / 255.0); cv::cvtColor(image, image, cv::COLOR_RGB2YUV); convert_mat(conv, image, denoise_level, scale, dst_w, dst_h, block_size, w2xc::IMAGE_Y); cv::cvtColor(image, image, cv::COLOR_YUV2RGB); image.convertTo(dsti, CV_8U, 255.0); } return 0; }
int w2xconv_convert_yuv(struct W2XConv *conv, unsigned char *dst, size_t dst_step_byte, /* float32x3 (src_w*ratio, src_h*ratio) */ unsigned char *src, size_t src_step_byte, /* float32x3 (src_w, src_h) */ int src_w, int src_h, int denoise_level, /* 0:none, 1:L1 denoise, other:L2 denoise */ double scale, int block_size) { int dst_h = src_h * scale; int dst_w = src_w * scale; bool is_rgb = (conv->impl->scale2_models[0]->getNInputPlanes() == 3); if (is_rgb) { setError(conv, W2XCONV_ERROR_RGB_MODEL_MISMATCH_TO_Y); return -1; } cv::Mat srci(src_h, src_w, CV_32FC3, src, src_step_byte); cv::Mat dsti(dst_h, dst_w, CV_32FC3, dst, dst_step_byte); cv::Mat image = srci.clone(); convert_mat(conv, image, denoise_level, scale, dst_w, dst_h, block_size, w2xc::IMAGE_Y); image.copyTo(dsti); return 0; }
int gmx_cluster(int argc,char *argv[]) { static char *desc[] = { "g_cluster can cluster structures with several different methods.", "Distances between structures can be determined from a trajectory", "or read from an XPM matrix file with the [TT]-dm[tt] option.", "RMS deviation after fitting or RMS deviation of atom-pair distances", "can be used to define the distance between structures.[PAR]", "single linkage: add a structure to a cluster when its distance to any", "element of the cluster is less than [TT]cutoff[tt].[PAR]", "Jarvis Patrick: add a structure to a cluster when this structure", "and a structure in the cluster have each other as neighbors and", "they have a least [TT]P[tt] neighbors in common. The neighbors", "of a structure are the M closest structures or all structures within", "[TT]cutoff[tt].[PAR]", "Monte Carlo: reorder the RMSD matrix using Monte Carlo.[PAR]", "diagonalization: diagonalize the RMSD matrix.[PAR]" "gromos: use algorithm as described in Daura [IT]et al.[it]", "([IT]Angew. Chem. Int. Ed.[it] [BB]1999[bb], [IT]38[it], pp 236-240).", "Count number of neighbors using cut-off, take structure with", "largest number of neighbors with all its neighbors as cluster", "and eleminate it from the pool of clusters. Repeat for remaining", "structures in pool.[PAR]", "When the clustering algorithm assigns each structure to exactly one", "cluster (single linkage, Jarvis Patrick and gromos) and a trajectory", "file is supplied, the structure with", "the smallest average distance to the others or the average structure", "or all structures for each cluster will be written to a trajectory", "file. When writing all structures, separate numbered files are made", "for each cluster.[PAR]" "Two output files are always written:[BR]", "[TT]-o[tt] writes the RMSD values in the upper left half of the matrix", "and a graphical depiction of the clusters in the lower right half", "When [TT]-minstruct[tt] = 1 the graphical depiction is black", "when two structures are in the same cluster.", "When [TT]-minstruct[tt] > 1 different colors will be used for each", "cluster.[BR]", "[TT]-g[tt] writes information on the options used and a detailed list", "of all clusters and their members.[PAR]", "Additionally, a number of optional output files can be written:[BR]", "[TT]-dist[tt] writes the RMSD distribution.[BR]", "[TT]-ev[tt] writes the eigenvectors of the RMSD matrix", "diagonalization.[BR]", "[TT]-sz[tt] writes the cluster sizes.[BR]", "[TT]-tr[tt] writes a matrix of the number transitions between", "cluster pairs.[BR]", "[TT]-ntr[tt] writes the total number of transitions to or from", "each cluster.[BR]", "[TT]-clid[tt] writes the cluster number as a function of time.[BR]", "[TT]-cl[tt] writes average (with option [TT]-av[tt]) or central", "structure of each cluster or writes numbered files with cluster members", "for a selected set of clusters (with option [TT]-wcl[tt], depends on", "[TT]-nst[tt] and [TT]-rmsmin[tt]).[BR]", }; FILE *fp,*log; int i,i1,i2,j,nf,nrms; matrix box; rvec *xtps,*usextps,*x1,**xx=NULL; char *fn,*trx_out_fn; t_clusters clust; t_mat *rms; real *eigval; t_topology top; int ePBC; t_atoms useatoms; t_matrix *readmat; real *tmp; int isize=0,ifsize=0,iosize=0; atom_id *index=NULL, *fitidx, *outidx; char *grpname; real rmsd,**d1,**d2,*time,time_invfac,*mass=NULL; char buf[STRLEN],buf1[80],title[STRLEN]; bool bAnalyze,bUseRmsdCut,bJP_RMSD=FALSE,bReadMat,bReadTraj; int method,ncluster=0; static char *methodname[] = { NULL, "linkage", "jarvis-patrick","monte-carlo", "diagonalization", "gromos", NULL }; enum { m_null, m_linkage, m_jarvis_patrick, m_monte_carlo, m_diagonalize, m_gromos, m_nr }; /* Set colors for plotting: white = zero RMS, black = maximum */ static t_rgb rlo_top = { 1.0, 1.0, 1.0 }; static t_rgb rhi_top = { 0.0, 0.0, 0.0 }; static t_rgb rlo_bot = { 1.0, 1.0, 1.0 }; static t_rgb rhi_bot = { 0.0, 0.0, 1.0 }; static int nlevels=40,skip=1; static real scalemax=-1.0,rmsdcut=0.1,rmsmin=0.0; static bool bRMSdist=FALSE,bBinary=FALSE,bAverage=FALSE,bFit=TRUE; static int niter=10000,seed=1993,write_ncl=0,write_nst=1,minstruct=1; static real kT=1e-3; static int M=10,P=3; t_pargs pa[] = { { "-dista", FALSE, etBOOL, {&bRMSdist}, "Use RMSD of distances instead of RMS deviation" }, { "-nlevels",FALSE,etINT, {&nlevels}, "Discretize RMSD matrix in # levels" }, { "-cutoff",FALSE, etREAL, {&rmsdcut}, "RMSD cut-off (nm) for two structures to be neighbor" }, { "-fit", FALSE, etBOOL, {&bFit}, "Use least squares fitting before RMSD calculation" }, { "-max", FALSE, etREAL, {&scalemax}, "Maximum level in RMSD matrix" }, { "-skip", FALSE, etINT, {&skip}, "Only analyze every nr-th frame" }, { "-av", FALSE, etBOOL, {&bAverage}, "Write average iso middle structure for each cluster" }, { "-wcl", FALSE, etINT, {&write_ncl}, "Write all structures for first # clusters to numbered files" }, { "-nst", FALSE, etINT, {&write_nst}, "Only write all structures if more than # per cluster" }, { "-rmsmin",FALSE, etREAL, {&rmsmin}, "minimum rms difference with rest of cluster for writing structures" }, { "-method",FALSE, etENUM, {methodname}, "Method for cluster determination" }, { "-minstruct", FALSE, etINT, {&minstruct}, "Minimum number of structures in cluster for coloring in the xpm file" }, { "-binary",FALSE, etBOOL, {&bBinary}, "Treat the RMSD matrix as consisting of 0 and 1, where the cut-off " "is given by -cutoff" }, { "-M", FALSE, etINT, {&M}, "Number of nearest neighbors considered for Jarvis-Patrick algorithm, " "0 is use cutoff" }, { "-P", FALSE, etINT, {&P}, "Number of identical nearest neighbors required to form a cluster" }, { "-seed", FALSE, etINT, {&seed}, "Random number seed for Monte Carlo clustering algorithm" }, { "-niter", FALSE, etINT, {&niter}, "Number of iterations for MC" }, { "-kT", FALSE, etREAL, {&kT}, "Boltzmann weighting factor for Monte Carlo optimization " "(zero turns off uphill steps)" } }; t_filenm fnm[] = { { efTRX, "-f", NULL, ffOPTRD }, { efTPS, "-s", NULL, ffOPTRD }, { efNDX, NULL, NULL, ffOPTRD }, { efXPM, "-dm", "rmsd", ffOPTRD }, { efXPM, "-o", "rmsd-clust", ffWRITE }, { efLOG, "-g", "cluster", ffWRITE }, { efXVG, "-dist", "rmsd-dist", ffOPTWR }, { efXVG, "-ev", "rmsd-eig", ffOPTWR }, { efXVG, "-sz", "clust-size", ffOPTWR}, { efXPM, "-tr", "clust-trans",ffOPTWR}, { efXVG, "-ntr", "clust-trans",ffOPTWR}, { efXVG, "-clid", "clust-id.xvg",ffOPTWR}, { efTRX, "-cl", "clusters.pdb", ffOPTWR } }; #define NFILE asize(fnm) CopyRight(stderr,argv[0]); parse_common_args(&argc,argv,PCA_CAN_VIEW | PCA_CAN_TIME | PCA_TIME_UNIT | PCA_BE_NICE, NFILE,fnm,asize(pa),pa,asize(desc),desc,0,NULL); /* parse options */ bReadMat = opt2bSet("-dm",NFILE,fnm); bReadTraj = opt2bSet("-f",NFILE,fnm) || !bReadMat; if ( opt2parg_bSet("-av",asize(pa),pa) || opt2parg_bSet("-wcl",asize(pa),pa) || opt2parg_bSet("-nst",asize(pa),pa) || opt2parg_bSet("-rmsmin",asize(pa),pa) || opt2bSet("-cl",NFILE,fnm) ) trx_out_fn = opt2fn("-cl",NFILE,fnm); else trx_out_fn = NULL; if (bReadMat && time_factor()!=1) { fprintf(stderr, "\nWarning: assuming the time unit in %s is %s\n", opt2fn("-dm",NFILE,fnm),time_unit()); } if (trx_out_fn && !bReadTraj) fprintf(stderr,"\nWarning: " "cannot write cluster structures without reading trajectory\n" " ignoring option -cl %s\n", trx_out_fn); method=1; while ( method < m_nr && strcasecmp(methodname[0], methodname[method])!=0 ) method++; if (method == m_nr) gmx_fatal(FARGS,"Invalid method"); bAnalyze = (method == m_linkage || method == m_jarvis_patrick || method == m_gromos ); /* Open log file */ log = ftp2FILE(efLOG,NFILE,fnm,"w"); fprintf(stderr,"Using %s method for clustering\n",methodname[0]); fprintf(log,"Using %s method for clustering\n",methodname[0]); /* check input and write parameters to log file */ bUseRmsdCut = FALSE; if (method == m_jarvis_patrick) { bJP_RMSD = (M == 0) || opt2parg_bSet("-cutoff",asize(pa),pa); if ((M<0) || (M == 1)) gmx_fatal(FARGS,"M (%d) must be 0 or larger than 1",M); if (M < 2) { sprintf(buf1,"Will use P=%d and RMSD cutoff (%g)",P,rmsdcut); bUseRmsdCut = TRUE; } else { if (P >= M) gmx_fatal(FARGS,"Number of neighbors required (P) must be less than M"); if (bJP_RMSD) { sprintf(buf1,"Will use P=%d, M=%d and RMSD cutoff (%g)",P,M,rmsdcut); bUseRmsdCut = TRUE; } else sprintf(buf1,"Will use P=%d, M=%d",P,M); } ffprintf1(stderr,log,buf,"%s for determining the neighbors\n\n",buf1); } else /* method != m_jarvis */ bUseRmsdCut = ( bBinary || method == m_linkage || method == m_gromos ); if (bUseRmsdCut && method != m_jarvis_patrick) fprintf(log,"Using RMSD cutoff %g nm\n",rmsdcut); if ( method==m_monte_carlo ) fprintf(log,"Using %d iterations\n",niter); if (skip < 1) gmx_fatal(FARGS,"skip (%d) should be >= 1",skip); /* get input */ if (bReadTraj) { /* don't read mass-database as masses (and top) are not used */ read_tps_conf(ftp2fn(efTPS,NFILE,fnm),buf,&top,&ePBC,&xtps,NULL,box, bAnalyze); fprintf(stderr,"\nSelect group for least squares fit%s:\n", bReadMat?"":" and RMSD calculation"); get_index(&(top.atoms),ftp2fn_null(efNDX,NFILE,fnm), 1,&ifsize,&fitidx,&grpname); if (trx_out_fn) { fprintf(stderr,"\nSelect group for output:\n"); get_index(&(top.atoms),ftp2fn_null(efNDX,NFILE,fnm), 1,&iosize,&outidx,&grpname); /* merge and convert both index groups: */ /* first copy outidx to index. let outidx refer to elements in index */ snew(index,iosize); isize = iosize; for(i=0; i<iosize; i++) { index[i]=outidx[i]; outidx[i]=i; } /* now lookup elements from fitidx in index, add them if necessary and also let fitidx refer to elements in index */ for(i=0; i<ifsize; i++) { j=0; while (j<isize && index[j]!=fitidx[i]) j++; if (j>=isize) { /* slow this way, but doesn't matter much */ isize++; srenew(index,isize); } index[j]=fitidx[i]; fitidx[i]=j; } } else { /* !trx_out_fn */ isize = ifsize; snew(index, isize); for(i=0; i<ifsize; i++) { index[i]=fitidx[i]; fitidx[i]=i; } } } /* Initiate arrays */ snew(d1,isize); snew(d2,isize); for(i=0; (i<isize); i++) { snew(d1[i],isize); snew(d2[i],isize); } if (bReadTraj) { /* Loop over first coordinate file */ fn = opt2fn("-f",NFILE,fnm); xx = read_whole_trj(fn,isize,index,skip,&nf,&time); convert_times(nf, time); if (!bRMSdist || bAnalyze) { /* Center all frames on zero */ snew(mass,isize); for(i=0; i<ifsize; i++) mass[fitidx[i]] = top.atoms.atom[index[fitidx[i]]].m; if (bFit) for(i=0; i<nf; i++) reset_x(ifsize,fitidx,isize,NULL,xx[i],mass); } } if (bReadMat) { fprintf(stderr,"Reading rms distance matrix "); read_xpm_matrix(opt2fn("-dm",NFILE,fnm),&readmat); fprintf(stderr,"\n"); if (readmat[0].nx != readmat[0].ny) gmx_fatal(FARGS,"Matrix (%dx%d) is not square", readmat[0].nx,readmat[0].ny); if (bReadTraj && bAnalyze && (readmat[0].nx != nf)) gmx_fatal(FARGS,"Matrix size (%dx%d) does not match the number of " "frames (%d)",readmat[0].nx,readmat[0].ny,nf); nf = readmat[0].nx; sfree(time); time = readmat[0].axis_x; time_invfac = time_invfactor(); for(i=0; i<nf; i++) time[i] *= time_invfac; rms = init_mat(readmat[0].nx,method == m_diagonalize); convert_mat(&(readmat[0]),rms); nlevels = readmat[0].nmap; } else { /* !bReadMat */ rms = init_mat(nf,method == m_diagonalize); nrms = (nf*(nf-1))/2; if (!bRMSdist) { fprintf(stderr,"Computing %dx%d RMS deviation matrix\n",nf,nf); snew(x1,isize); for(i1=0; (i1<nf); i1++) { for(i2=i1+1; (i2<nf); i2++) { for(i=0; i<isize; i++) copy_rvec(xx[i1][i],x1[i]); if (bFit) do_fit(isize,mass,xx[i2],x1); rmsd = rmsdev(isize,mass,xx[i2],x1); set_mat_entry(rms,i1,i2,rmsd); } nrms -= (nf-i1-1); fprintf(stderr,"\r# RMSD calculations left: %d ",nrms); } } else { /* bRMSdist */ fprintf(stderr,"Computing %dx%d RMS distance deviation matrix\n",nf,nf); for(i1=0; (i1<nf); i1++) { calc_dist(isize,xx[i1],d1); for(i2=i1+1; (i2<nf); i2++) { calc_dist(isize,xx[i2],d2); set_mat_entry(rms,i1,i2,rms_dist(isize,d1,d2)); } nrms -= (nf-i1-1); fprintf(stderr,"\r# RMSD calculations left: %d ",nrms); } } fprintf(stderr,"\n\n"); } ffprintf2(stderr,log,buf,"The RMSD ranges from %g to %g nm\n", rms->minrms,rms->maxrms); ffprintf1(stderr,log,buf,"Average RMSD is %g\n",2*rms->sumrms/(nf*(nf-1))); ffprintf1(stderr,log,buf,"Number of structures for matrix %d\n",nf); ffprintf1(stderr,log,buf,"Energy of the matrix is %g nm\n",mat_energy(rms)); if (bUseRmsdCut && (rmsdcut < rms->minrms || rmsdcut > rms->maxrms) ) fprintf(stderr,"WARNING: rmsd cutoff %g is outside range of rmsd values " "%g to %g\n",rmsdcut,rms->minrms,rms->maxrms); if (bAnalyze && (rmsmin < rms->minrms) ) fprintf(stderr,"WARNING: rmsd minimum %g is below lowest rmsd value %g\n", rmsmin,rms->minrms); if (bAnalyze && (rmsmin > rmsdcut) ) fprintf(stderr,"WARNING: rmsd minimum %g is above rmsd cutoff %g\n", rmsmin,rmsdcut); /* Plot the rmsd distribution */ rmsd_distribution(opt2fn("-dist",NFILE,fnm),rms); if (bBinary) { for(i1=0; (i1 < nf); i1++) for(i2=0; (i2 < nf); i2++) if (rms->mat[i1][i2] < rmsdcut) rms->mat[i1][i2] = 0; else rms->mat[i1][i2] = 1; } snew(clust.cl,nf); switch (method) { case m_linkage: /* Now sort the matrix and write it out again */ gather(rms,rmsdcut,&clust); break; case m_diagonalize: /* Do a diagonalization */ snew(eigval,nf); snew(tmp,nf*nf); memcpy(tmp,rms->mat[0],nf*nf*sizeof(real)); eigensolver(tmp,nf,0,nf,eigval,rms->mat[0]); sfree(tmp); fp = xvgropen(opt2fn("-ev",NFILE,fnm),"RMSD matrix Eigenvalues", "Eigenvector index","Eigenvalues (nm\\S2\\N)"); for(i=0; (i<nf); i++) fprintf(fp,"%10d %10g\n",i,eigval[i]); ffclose(fp); break; case m_monte_carlo: mc_optimize(log,rms,niter,&seed,kT); swap_mat(rms); reset_index(rms); break; case m_jarvis_patrick: jarvis_patrick(rms->nn,rms->mat,M,P,bJP_RMSD ? rmsdcut : -1,&clust); break; case m_gromos: gromos(rms->nn,rms->mat,rmsdcut,&clust); break; default: gmx_fatal(FARGS,"DEATH HORROR unknown method \"%s\"",methodname[0]); } if (method == m_monte_carlo || method == m_diagonalize) fprintf(stderr,"Energy of the matrix after clustering is %g nm\n", mat_energy(rms)); if (bAnalyze) { if (minstruct > 1) { ncluster = plot_clusters(nf,rms->mat,&clust,nlevels,minstruct); } else { mark_clusters(nf,rms->mat,rms->maxrms,&clust); } init_t_atoms(&useatoms,isize,FALSE); snew(usextps, isize); useatoms.resname=top.atoms.resname; for(i=0; i<isize; i++) { useatoms.atomname[i]=top.atoms.atomname[index[i]]; useatoms.atom[i].resnr=top.atoms.atom[index[i]].resnr; useatoms.nres=max(useatoms.nres,useatoms.atom[i].resnr+1); copy_rvec(xtps[index[i]],usextps[i]); } useatoms.nr=isize; analyze_clusters(nf,&clust,rms->mat,isize,&useatoms,usextps,mass,xx,time, ifsize,fitidx,iosize,outidx, bReadTraj?trx_out_fn:NULL, opt2fn_null("-sz",NFILE,fnm), opt2fn_null("-tr",NFILE,fnm), opt2fn_null("-ntr",NFILE,fnm), opt2fn_null("-clid",NFILE,fnm), bAverage, write_ncl, write_nst, rmsmin, bFit, log, rlo_bot,rhi_bot); } ffclose(log); if (bBinary && !bAnalyze) /* Make the clustering visible */ for(i2=0; (i2 < nf); i2++) for(i1=i2+1; (i1 < nf); i1++) if (rms->mat[i1][i2]) rms->mat[i1][i2] = rms->maxrms; fp = opt2FILE("-o",NFILE,fnm,"w"); fprintf(stderr,"Writing rms distance/clustering matrix "); if (bReadMat) { write_xpm(fp,0,readmat[0].title,readmat[0].legend,readmat[0].label_x, readmat[0].label_y,nf,nf,readmat[0].axis_x,readmat[0].axis_y, rms->mat,0.0,rms->maxrms,rlo_top,rhi_top,&nlevels); } else { sprintf(buf,"Time (%s)",time_unit()); sprintf(title,"RMS%sDeviation / Cluster Index", bRMSdist ? " Distance " : " "); if (minstruct > 1) { write_xpm_split(fp,0,title,"RMSD (nm)",buf,buf, nf,nf,time,time,rms->mat,0.0,rms->maxrms,&nlevels, rlo_top,rhi_top,0.0,(real) ncluster, &ncluster,TRUE,rlo_bot,rhi_bot); } else { write_xpm(fp,0,title,"RMSD (nm)",buf,buf, nf,nf,time,time,rms->mat,0.0,rms->maxrms, rlo_top,rhi_top,&nlevels); } } fprintf(stderr,"\n"); ffclose(fp); /* now show what we've done */ do_view(opt2fn("-o",NFILE,fnm),"-nxy"); do_view(opt2fn_null("-sz",NFILE,fnm),"-nxy"); if (method == m_diagonalize) do_view(opt2fn_null("-ev",NFILE,fnm),"-nxy"); do_view(opt2fn("-dist",NFILE,fnm),"-nxy"); if (bAnalyze) { do_view(opt2fn_null("-tr",NFILE,fnm),"-nxy"); do_view(opt2fn_null("-ntr",NFILE,fnm),"-nxy"); do_view(opt2fn_null("-clid",NFILE,fnm),"-nxy"); } /* Thank the user for her patience */ thanx(stderr); return 0; }