/** For each point p, calculate function Kdist(p) which is the distance of * the Kth nearest point to p. */ void Cluster_DBSCAN::ComputeKdist( int Kval, std::vector<int> const& FramesToCluster ) const { std::vector<double> dists; std::vector<double> Kdist; dists.reserve( FramesToCluster.size() ); Kdist.reserve( FramesToCluster.size() ); std::string outfilename = k_prefix_ + "Kdist." + integerToString(Kval) + ".dat"; mprintf("\tDBSCAN: Calculating Kdist(%i), output to %s\n", Kval, outfilename.c_str()); for (std::vector<int>::const_iterator point = FramesToCluster.begin(); point != FramesToCluster.end(); ++point) { // Store distances from this point dists.clear(); for (std::vector<int>::const_iterator otherpoint = FramesToCluster.begin(); otherpoint != FramesToCluster.end(); ++otherpoint) dists.push_back( FrameDistances_.GetFdist(*point, *otherpoint) ); // Sort distances - first dist should always be 0 std::sort(dists.begin(), dists.end()); Kdist.push_back( dists[Kval] ); } std::sort( Kdist.begin(), Kdist.end() ); CpptrajFile Outfile; Outfile.OpenWrite(outfilename); Outfile.Printf("%-8s %1i%-11s\n", "#Point", Kval,"-dist"); // Write out largest to smallest unsigned int ik = 0; for (std::vector<double>::reverse_iterator k = Kdist.rbegin(); k != Kdist.rend(); ++k, ++ik) Outfile.Printf("%8u %12.4f\n", ik, *k); Outfile.CloseFile(); }
int DataIO_OpenDx::WriteGrid(DataSet const& setIn, CpptrajFile& outfile) const { DataSet_3D const& set = static_cast<DataSet_3D const&>( setIn ); Vec3 oxyz = set.Bin().GridOrigin(); if (gridWriteMode_ == BIN_CENTER) // Origin needs to be shifted to center of bin located at 0,0,0 oxyz = set.Bin().Center(0,0,0); // Print the OpenDX header WriteDxHeader(outfile, set.NX(), set.NY(), set.NZ(), set.NX(), set.NY(), set.NZ(), set.Bin().Ucell(), oxyz); // Now print out the data. size_t gridsize = set.Size(); if (gridsize == 1) outfile.Printf("%g\n", set[0]); else if (gridsize == 2) outfile.Printf("%g %g\n", set[0], set[1]); else if (gridsize > 2) { // Data is already in row-major form (z-axis changes // fastest), so no need to do any kind of data adjustment for (size_t i = 0UL; i < gridsize - 2UL; i += 3UL) outfile.Printf("%g %g %g\n", set[i], set[i+1], set[i+2]); // Print out any points we may have missed switch (gridsize % 3) { case 2: outfile.Printf("%g %g\n", set[gridsize-2], set[gridsize-1]); break; case 1: outfile.Printf("%g\n", set[gridsize-1]); break; } } return 0; }
void DataSet_GridDbl::WriteBuffer(CpptrajFile& outfile, SizeArray const& pIn) const { size_t x = pIn[0]; size_t y = pIn[1]; size_t z = pIn[2]; if ( x >= grid_.NX() || y >= grid_.NY() || z >= grid_.NZ() ) outfile.Printf(format_.fmt(), 0.0); else outfile.Printf(format_.fmt(), grid_.element(x,y,z)); }
// DataSet_Vector::WriteBuffer() void DataSet_Vector::WriteBuffer(CpptrajFile &cbuffer, SizeArray const& pIn) const { if (pIn[0] >= vectors_.size()) { cbuffer.Printf(format_.fmt(), 0.0, 0.0, 0.0, 0.0, 0.0, 0.0); // VXYZ OXYZ } else { Vec3 const& Vxyz = vectors_[pIn[0]]; Vec3 const& Oxyz = OXYZ(pIn[0]); cbuffer.Printf(format_.fmt(), Vxyz[0], Vxyz[1], Vxyz[2], Oxyz[0], Oxyz[1], Oxyz[2]); } }
void TriangleMatrix::Write2D( CpptrajFile& outfile, int xIn, int yIn ) { size_t x = (size_t)xIn; size_t y = (size_t)yIn; if ( xIn==yIn || xIn < 0 || yIn < 0 || x >= nrows_ || y >= nrows_ ) outfile.Printf(data_format_, 0.0); else { size_t index = calcIndex(x, y); outfile.Printf(data_format_, elements_[index]); } }
void Cluster_DPeaks::ClusterResults(CpptrajFile& outfile) const { outfile.Printf("#Algorithm: DPeaks epsilon %g\n", epsilon_); if (calc_noise_) { outfile.Printf("#NOISE_FRAMES:"); std::vector<int> noiseFrames; for (Carray::const_iterator point = Points_.begin(); point != Points_.end(); ++point) if (point->Cnum() == -1) noiseFrames.push_back( point->Fnum()+1 ); std::sort( noiseFrames.begin(), noiseFrames.end() ); for (std::vector<int>::const_iterator f = noiseFrames.begin(); f != noiseFrames.end(); ++f) outfile.Printf(" %i", *f); outfile.Printf("\n"); outfile.Printf("#Number_of_noise_frames: %zu\n", noiseFrames.size()); } }
// WriteNameToBuffer() void DataIO_Std::WriteNameToBuffer(CpptrajFile& fileIn, std::string const& label, int width, bool isLeftCol) { std::string temp_name = label; // If left aligning, add '#' to name; if (isLeftCol) { if (temp_name[0]!='#') { temp_name.insert(0,"#"); // Ensure that name will not be larger than column width. if ((int)temp_name.size() > width) temp_name.resize( width ); } } // Replace any spaces with underscores for (std::string::iterator tc = temp_name.begin(); tc != temp_name.end(); ++tc) if ( *tc == ' ' ) *tc = '_'; if (width >= (int)CpptrajFile::BUF_SIZE) // Protect against CpptrajFile buffer overflow fileIn.Write(temp_name.c_str(), temp_name.size()); else { // Set up header format string TextFormat::AlignType align; if (isLeftCol) align = TextFormat::LEFT; else align = TextFormat::RIGHT; TextFormat header_format(TextFormat::STRING, width, align); fileIn.Printf(header_format.fmt(), temp_name.c_str()); } }
// DataIO_Std::WriteByGroup() int DataIO_Std::WriteByGroup(CpptrajFile& file, DataSetList const& SetList, GroupType gtype) { int err = 0; bool firstWrite = true; DataSetList tmpdsl; std::vector<bool> setIsWritten(SetList.size(), false); unsigned int startIdx = 0; unsigned int nWritten = 0; while (nWritten < SetList.size()) { std::string currentName; Dimension currentDim; int currentNum = -1; switch (gtype) { case BY_NAME : currentName = SetList[startIdx]->Meta().Name(); break; case BY_ASPECT : currentName = SetList[startIdx]->Meta().Aspect(); break; case BY_IDX : currentNum = SetList[startIdx]->Meta().Idx(); break; case BY_ENS : currentNum = SetList[startIdx]->Meta().EnsembleNum(); break; case BY_DIM : currentDim = SetList[startIdx]->Dim(0); break; case NO_TYPE : return 1; } int firstNonMatch = -1; for (unsigned int idx = startIdx; idx != SetList.size(); idx++) { if (!setIsWritten[idx]) { bool match = false; switch (gtype) { case BY_NAME : match = (currentName == SetList[idx]->Meta().Name()); break; case BY_ASPECT : match = (currentName == SetList[idx]->Meta().Aspect()); break; case BY_IDX : match = (currentNum == SetList[idx]->Meta().Idx()); break; case BY_ENS : match = (currentNum == SetList[idx]->Meta().EnsembleNum()); break; case BY_DIM : match = (currentDim == SetList[idx]->Dim(0)); break; case NO_TYPE : return 1; } if (match) { tmpdsl.AddCopyOfSet( SetList[idx] ); setIsWritten[idx] = true; nWritten++; } else if (firstNonMatch == -1) firstNonMatch = (int)idx; } } if (firstNonMatch > -1) startIdx = (unsigned int)firstNonMatch; if (!firstWrite) file.Printf("\n"); else firstWrite = false; if (isInverted_) err += WriteDataInverted(file, tmpdsl); else err += WriteDataNormal(file, tmpdsl); tmpdsl.ClearAll(); } return err; }
// Cluster_DBSCAN::ClusterResults() void Cluster_DBSCAN::ClusterResults(CpptrajFile& outfile) const { outfile.Printf("#Algorithm: DBSCAN minpoints %i epsilon %g sieveToCentroid %i\n", minPoints_, epsilon_, (int)sieveToCentroid_); // List the number of noise points. outfile.Printf("#NOISE_FRAMES:"); unsigned int frame = 1; unsigned int numNoise = 0; for (std::vector<char>::const_iterator stat = Status_.begin(); stat != Status_.end(); ++stat, ++frame) { if ( *stat == NOISE ) { outfile.Printf(" %i", frame); ++numNoise; } } outfile.Printf("\n"); outfile.Printf("#Number_of_noise_frames: %u\n", numNoise); }
// DataIO_Std::WriteData3D() int DataIO_Std::WriteData3D( CpptrajFile& file, DataSetList const& setList) { int err = 0; for (DataSetList::const_iterator set = setList.begin(); set != setList.end(); ++set) { if (set != setList.begin()) file.Printf("\n"); err += WriteSet3D( *(*set), file ); } return err; }
int DataIO_OpenDx::WriteGridWrap(DataSet const& setIn, CpptrajFile& outfile) const { DataSet_3D const& set = static_cast<DataSet_3D const&>( setIn ); // Need to construct a grid mesh around bins, with points centered on the bins. int mesh_x = set.NX(); int mesh_y = set.NY(); int mesh_z = set.NZ(); // Origin needs to be shifted half grid spacing, i.e. it is the center of the // bin located at -1, -1, -1. Vec3 oxyz = set.Bin().Center(-1, -1, -1); // Print the OpenDX header WriteDxHeader(outfile, mesh_x+2, mesh_y+2, mesh_z+2, mesh_x, mesh_y, mesh_z, set.Bin().Ucell(), oxyz); // Print out the data. Start at bin -1, end on bin N. int nvals = 0; // Keep track of how many values printed on current line. if (gridWriteMode_ == WRAP) { int bi, bj, bk; for (int ii = -1; ii <= mesh_x; ++ii) { if (ii < 0 ) bi = mesh_x - 1; else if (ii == mesh_x) bi = 0; else bi = ii; for (int ij = -1; ij <= mesh_y; ++ij) { if (ij < 0 ) bj = mesh_y - 1; else if (ij == mesh_y) bj = 0; else bj = ij; for (int ik = -1; ik <= mesh_z; ++ik) { if (ik < 0 ) bk = mesh_z - 1; else if (ik == mesh_z) bk = 0; else bk = ik; outfile.Printf(" %g", set.GetElement(bi, bj, bk)); ++nvals; if (nvals == 5) { outfile.Printf("\n"); nvals = 0; } } } } } else { // EXTENDED for (int ii = -1; ii <= mesh_x; ++ii) { bool zero_x = (ii < 0 || ii == mesh_x); for (int ij = -1; ij <= mesh_y; ++ij) { bool zero_y = (ij < 0 || ij == mesh_y); for (int ik = -1; ik <= mesh_z; ++ik) { if (zero_x || zero_y || ik < 0 || ik == mesh_z) outfile.Printf(" 0"); else outfile.Printf(" %g", set.GetElement(ii, ij, ik)); ++nvals; if (nvals == 5) { outfile.Printf("\n"); nvals = 0; } } } } } if (nvals > 0) outfile.Printf("\n"); return 0; }
void Action_Pairwise::Print() { if (nframes_ < 1) return; // Divide matrices by # of frames double norm = 1.0 / (double)nframes_; for (unsigned int i = 0; i != vdwMat_->Size(); i++) { (*vdwMat_)[i] *= norm; (*eleMat_)[i] *= norm; } // Write out final results CpptrajFile AvgOut; if (AvgOut.OpenWrite( avgout_ )) return; if (nb_calcType_ == NORMAL) mprintf(" PAIRWISE: Writing all pairs with |<evdw>| > %.4f, |<eelec>| > %.4f\n", cut_evdw_, cut_eelec_); else if (nb_calcType_ == COMPARE_REF) mprintf(" PAIRWISE: Writing all pairs with |<dEvdw>| > %.4f, |<dEelec>| > %.4f\n", cut_evdw_, cut_eelec_); AvgOut.Printf("%-16s %5s -- %16s %5s : ENE\n","#Name1", "At1", "Name2", "At2"); for (AtomMask::const_iterator m1 = Mask0_.begin(); m1 != Mask0_.end(); ++m1) { for (AtomMask::const_iterator m2 = m1 + 1; m2 != Mask0_.end(); ++m2) { double EV = vdwMat_->GetElement(*m1, *m2); double EE = eleMat_->GetElement(*m1, *m2); bool outputv = ( fabs(EV) > cut_evdw_ ); bool outpute = ( fabs(EE) > cut_eelec_ ); if (outputv || outpute) { AvgOut.Printf("%16s %5i -- %16s %5i :", CurrentParm_->TruncResAtomName(*m1).c_str(), *m1 + 1, CurrentParm_->TruncResAtomName(*m2).c_str(), *m2 + 1); if (outputv) AvgOut.Printf(" EVDW= %12.5e", EV); if (outpute) AvgOut.Printf(" EELEC= %12.5e", EE); AvgOut.Printf("\n"); } } } }
void DataIO_OpenDx::WriteDxHeader(CpptrajFile& outfile, size_t NX, size_t NY, size_t NZ, double LX, double LY, double LZ, Matrix_3x3 const& ucell, Vec3 const& oxyz) const { outfile.Printf("object 1 class gridpositions counts %zu %zu %zu\n" "origin %g %g %g\ndelta %g %g %g\ndelta %g %g %g\ndelta %g %g %g\n" "object 2 class gridconnections counts %zu %zu %zu\n" "object 3 class array type double rank 0 items %zu data follows\n", NX, NY, NZ, oxyz[0], oxyz[1], oxyz[2], ucell[0]/LX, ucell[1]/LX, ucell[2]/LX, ucell[3]/LY, ucell[4]/LY, ucell[5]/LY, ucell[6]/LZ, ucell[7]/LZ, ucell[8]/LZ, NX, NY, NZ, NX*NY*NZ); }
// DataIO_Std::WriteCmatrix() int DataIO_Std::WriteCmatrix(CpptrajFile& file, DataSetList const& Sets) { for (DataSetList::const_iterator ds = Sets.begin(); ds != Sets.end(); ++ds) { if ( (*ds)->Group() != DataSet::CLUSTERMATRIX) { mprinterr("Error: Write of cluster matrix and other sets to same file not supported.\n" "Error: Skipping '%s'\n", (*ds)->legend()); continue; } DataSet_Cmatrix const& cm = static_cast<DataSet_Cmatrix const&>( *(*ds) ); int nrows = cm.OriginalNframes(); int col_width = std::max(3, DigitWidth( nrows ) + 1); int dat_width = std::max(cm.Format().Width(), (int)cm.Meta().Legend().size()) + 1; WriteNameToBuffer(file, "F1", col_width, true); WriteNameToBuffer(file, "F2", col_width, false); WriteNameToBuffer(file, cm.Meta().Legend(), dat_width, false); if (cm.SieveType() != ClusterSieve::NONE) file.Printf(" nframes %i", cm.OriginalNframes()); file.Printf("\n"); TextFormat col_fmt(TextFormat::INTEGER, col_width); TextFormat dat_fmt = cm.Format(); dat_fmt.SetFormatAlign(TextFormat::RIGHT); dat_fmt.SetFormatWidth( dat_width ); std::string total_fmt = col_fmt.Fmt() + col_fmt.Fmt() + dat_fmt.Fmt() + "\n"; //mprintf("DEBUG: format '%s'\n", total_fmt.c_str()); ClusterSieve::SievedFrames const& frames = cm.FramesToCluster(); int ntotal = (int)frames.size(); for (int idx1 = 0; idx1 != ntotal; idx1++) { int row = frames[idx1]; for (int idx2 = idx1 + 1; idx2 != ntotal; idx2++) { int col = frames[idx2]; file.Printf(total_fmt.c_str(), row+1, col+1, cm.GetFdist(col, row)); } } } return 0; }
// DataIO_OpenDx::WriteSet3D() int DataIO_OpenDx::WriteSet3D(DataSet const& setIn, CpptrajFile& outfile) const { if (setIn.Ndim() != 3) { mprinterr("Internal Error: DataSet %s in DataFile %s has %zu dimensions, expected 3.\n", setIn.legend(), outfile.Filename().full(), setIn.Ndim()); return 1; } int err = 0; switch ( gridWriteMode_ ) { case BIN_CORNER: case BIN_CENTER: err = WriteGrid( setIn, outfile ); break; case WRAP: case EXTENDED : err = WriteGridWrap( setIn, outfile ); break; } // Print tail if (err == 0) { // TODO: Make this an option //if (mode_ == CENTER) // outfile.Printf("\nobject \"density (%s) [A^-3]\" class field\n", // centerMask_.MaskString()); //else outfile.Printf("\nobject \"density [A^-3]\" class field\n"); } return err; }
// DataIO_Std::WriteDataInverted() int DataIO_Std::WriteDataInverted(CpptrajFile& file, DataSetList const& Sets) { if (Sets.empty() || CheckAllDims(Sets, 1)) return 1; // Determine size of largest DataSet. size_t maxFrames = DetermineMax( Sets ); // Write each set to a line. DataSet::SizeArray positions(1); // Set up x column format DataSetList::const_iterator set = Sets.begin(); TextFormat x_col_format; if (hasXcolumn_) x_col_format = XcolFmt(); else x_col_format = (*set)->Format(); for (; set != Sets.end(); ++set) { // Write dataset name as first column. WriteNameToBuffer( file, (*set)->Meta().Legend(), x_col_format.ColumnWidth(), false); // Write each frame to subsequent columns for (positions[0] = 0; positions[0] < maxFrames; positions[0]++) (*set)->WriteBuffer(file, positions); file.Printf("\n"); } return 0; }
// Analysis_Wavelet::Analyze() Analysis::RetType Analysis_Wavelet::Analyze() { // Step 1 - Create a matrix that is #atoms rows by #frames - 1 cols, // where matrix(frame, atom) is the distance that atom has // travelled from the previous frame. // TODO: Implement this in Action_Matrix()? mprintf(" WAVELET:\n"); // First set up atom mask. if (coords_->Top().SetupIntegerMask( mask_ )) return Analysis::ERR; mask_.MaskInfo(); int natoms = mask_.Nselected(); int nframes = (int)coords_->Size(); if (natoms < 1 || nframes < 2) { mprinterr("Error: Not enough frames (%i) or atoms (%i) in '%s'\n", nframes, natoms, coords_->legend()); return Analysis::ERR; } Matrix<double> d_matrix; mprintf("\t%i frames, %i atoms, distance matrix will require %.2f MB\n", (double)d_matrix.sizeInBytes(nframes, natoms) / (1024.0*1024.0)); d_matrix.resize(nframes, natoms); // Get initial frame. Frame currentFrame, lastFrame; currentFrame.SetupFrameFromMask( mask_, coords_->Top().Atoms() ); lastFrame = currentFrame; coords_->GetFrame( 0, lastFrame, mask_ ); // Iterate over frames for (int frm = 1; frm != nframes; frm++) { coords_->GetFrame( frm, currentFrame, mask_ ); int idx = frm; // Position in distance matrix; start at column 'frame' for (int at = 0; at != natoms; at++, idx += nframes) // Distance of atom in currentFrame from its position in lastFrame. d_matrix[idx] = sqrt(DIST2_NoImage( currentFrame.XYZ(at), lastFrame.XYZ(at) )); //lastFrame = currentFrame; // TODO: Re-enable? } # ifdef DEBUG_WAVELET // DEBUG: Write matrix to file. CpptrajFile dmatrixOut; // DEBUG dmatrixOut.OpenWrite("dmatrix.dat"); Matrix<double>::iterator mval = d_matrix.begin(); for (int row = 0; row != natoms; row++) { for (int col = 0; col != nframes; col++) dmatrixOut.Printf("%g ", *(mval++)); dmatrixOut.Printf("\n"); } dmatrixOut.CloseFile(); # endif // Precompute some factors for calculating scaled wavelets. const double one_over_sqrt_N = 1.0 / sqrt(static_cast<double>( nframes )); std::vector<int> arrayK( nframes ); arrayK[0] = -1 * (nframes/2); for (int i = 1; i != nframes; i++) arrayK[i] = arrayK[i-1] + 1; # ifdef DEBUG_WAVELET mprintf("DEBUG: K:"); for (std::vector<int>::const_iterator kval = arrayK.begin(); kval != arrayK.end(); ++kval) mprintf(" %i", *kval); mprintf("\n"); # endif // Step 2 - Get FFT of wavelet for each scale. PubFFT pubfft; pubfft.SetupFFTforN( nframes ); mprintf("\tMemory required for scaled wavelet array: %.2f MB\n", (double)(2 * nframes * nb_ * sizeof(double)) / (1024 * 1024)); typedef std::vector<ComplexArray> WaveletArray; WaveletArray FFT_of_Scaled_Wavelets; FFT_of_Scaled_Wavelets.reserve( nb_ ); typedef std::vector<double> Darray; Darray scaleVector; scaleVector.reserve( nb_ ); Darray MIN( nb_ * 2 ); for (int iscale = 0; iscale != nb_; iscale++) { // Calculate and store scale factor. scaleVector.push_back( S0_ * pow(2.0, iscale * ds_) ); // Populate MIN array MIN[iscale ] = (0.00647*pow((correction_*scaleVector.back()),1.41344)+19.7527)*chival_; MIN[iscale+nb_] = correction_*scaleVector.back(); // Calculate scalved wavelet ComplexArray scaledWavelet; switch (wavelet_type_) { case W_MORLET: scaledWavelet = F_Morlet(arrayK, scaleVector.back()); break; case W_PAUL : scaledWavelet = F_Paul(arrayK, scaleVector.back()); break; case W_NONE : return Analysis::ERR; // Sanity check } # ifdef DEBUG_WAVELET PrintComplex("wavelet_before_fft", scaledWavelet); # endif // Perform FFT pubfft.Forward( scaledWavelet ); // Normalize scaledWavelet.Normalize( one_over_sqrt_N ); # ifdef DEBUG_WAVELET PrintComplex("wavelet_after_fft", scaledWavelet); # endif FFT_of_Scaled_Wavelets.push_back( scaledWavelet ); } # ifdef DEBUG_WAVELET mprintf("DEBUG: Scaling factors:"); for (Darray::const_iterator sval = scaleVector.begin(); sval != scaleVector.end(); ++sval) mprintf(" %g", *sval); mprintf("\n"); mprintf("DEBUG: MIN:"); for (int i = 0; i != nb_; i++) mprintf(" %g", MIN[i]); mprintf("\n"); # endif // Step 3 - For each atom, calculate the convolution of scaled wavelets // with rows (atom distance vs frame) via dot product of the // frequency domains, i.e. Fourier-transformed, followed by an // inverse FT. DataSet_MatrixFlt& OUT = static_cast<DataSet_MatrixFlt&>( *output_ ); mprintf("\tMemory required for output matrix: %.2f MB\n", (double)Matrix<float>::sizeInBytes(nframes, natoms)/(1024.0*1024.0)); OUT.Allocate2D( nframes, natoms ); // Should initialize to zero Matrix<double> MAX; mprintf("\tMemory required for Max array: %.2f MB\n", (double)MAX.sizeInBytes(nframes, natoms)/(1024.0*1024.0)); MAX.resize( nframes, natoms ); Darray magnitude( nframes ); // Scratch space for calculating magnitude across rows for (int at = 0; at != natoms; at++) { ComplexArray AtomSignal( nframes ); // Initializes to zero // Calculate the distance variance for this atom and populate the array. int midx = at * nframes; // Index into d_matrix int cidx = 0; // Index into AtomSignal double d_avg = 0.0; double d_var = 0.0; for (int frm = 0; frm != nframes; frm++, cidx += 2, midx++) { d_avg += d_matrix[midx]; d_var += (d_matrix[midx] * d_matrix[midx]); AtomSignal[cidx] = d_matrix[midx]; } d_var = (d_var - ((d_avg * d_avg) / (double)nframes)) / ((double)(nframes - 1)); # ifdef DEBUG_WAVELET mprintf("VARIANCE: %g\n", d_var); # endif double var_norm = 1.0 / d_var; // Calculate FT of atom signal pubfft.Forward( AtomSignal ); # ifdef DEBUG_WAVELET PrintComplex("AtomSignal", AtomSignal); # endif // Normalize AtomSignal.Normalize( one_over_sqrt_N ); // Calculate dot product of atom signal with each scaled FT wavelet for (int iscale = 0; iscale != nb_; iscale++) { ComplexArray dot = AtomSignal.TimesComplexConj( FFT_of_Scaled_Wavelets[iscale] ); // Inverse FT of dot product pubfft.Back( dot ); # ifdef DEBUG_WAVELET PrintComplex("InverseFT_Dot", dot); # endif // Chi-squared testing midx = at * nframes; cidx = 0; for (int frm = 0; frm != nframes; frm++, cidx += 2, midx++) { magnitude[frm] = (dot[cidx]*dot[cidx] + dot[cidx+1]*dot[cidx+1]) * var_norm; if (magnitude[frm] < MIN[iscale]) magnitude[frm] = 0.0; if (magnitude[frm] > MAX[midx]) { MAX[midx] = magnitude[frm]; //Indices[midx] = iscale OUT[midx] = (float)(correction_ * scaleVector[iscale]); } } # ifdef DEBUG_WAVELET mprintf("DEBUG: AbsoluteValue:"); for (Darray::const_iterator dval = magnitude.begin(); dval != magnitude.end(); ++dval) mprintf(" %g", *dval); mprintf("\n"); # endif } // END loop over scales } // END loop over atoms # ifdef DEBUG_WAVELET // DEBUG: Print MAX CpptrajFile maxmatrixOut; // DEBUG maxmatrixOut.OpenWrite("maxmatrix.dat"); for (int col = 0; col != nframes; col++) { for (int row = 0; row != natoms; row++) maxmatrixOut.Printf("%g ", MAX.element(col, row)); maxmatrixOut.Printf("\n"); } maxmatrixOut.CloseFile(); # endif return Analysis::OK; }
int Cluster_DPeaks::ChoosePointsAutomatically() { // Right now all density values are discrete. Try to choose outliers at each // value for which there is density.; /* // For each point, calculate average distance (X,Y) to points in next and // previous density values. const double dens_cut = 3.0 * 3.0; const double dist_cut = 1.32 * 1.32; for (Carray::const_iterator point0 = Points_.begin(); point0 != Points_.end(); ++point0) { int Npts = 0; for (Carray::const_iterator point1 = Points_.begin(); point1 != Points_.end(); ++point1) { if (point0 != point1) { // Only do this for close densities double dX = (double)(point0->PointsWithinEps() - point1->PointsWithinEps()); double dX2 = dX * dX; double dY = (point0->Dist() - point1->Dist()); double dY2 = dY * dY; if (dX2 < dens_cut && dY2 < dist_cut) { Npts++; } } } mprintf("%i %i %i\n", point0->PointsWithinEps(), point0->Fnum()+1, Npts); } */ /* CpptrajFile tempOut; tempOut.OpenWrite("temp.dat"); int currentDensity = -1; double distAv = 0.0; double distSD = 0.0; double sumWts = 0.0; int nValues = 0; Carray::const_iterator lastPoint = Points_.end() + 1; for (Carray::const_iterator point = Points_.begin(); point != lastPoint; ++point) { if (point == Points_.end() || point->PointsWithinEps() != currentDensity) { if (nValues > 0) { distAv = distAv / sumWts; //(double)nValues; distSD = (distSD / sumWts) - (distAv * distAv); if (distSD > 0.0) distSD = sqrt(distSD); else distSD = 0.0; //mprintf("Density %i: %i values Avg= %g SD= %g SumWts= %g\n", currentDensity, // nValues, distAv, distSD, sumWts); tempOut.Printf("%i %g\n", currentDensity, distAv); } if (point == Points_.end()) break; currentDensity = point->PointsWithinEps(); distAv = 0.0; distSD = 0.0; sumWts = 0.0; nValues = 0; } double wt = exp(point->Dist()); double dval = point->Dist() * wt; sumWts += wt; distAv += dval; distSD += (dval * dval); nValues++; } tempOut.CloseFile(); */ // BEGIN CALCULATING WEIGHTED DISTANCE AVERAGE CpptrajFile tempOut; tempOut.OpenWrite("temp.dat"); DataSet_Mesh weightedAverage; Carray::const_iterator cp = Points_.begin(); // Skip local density of 0. //while (cp->PointsWithinEps() == 0 && cp != Points_.end()) ++cp; while (cp != Points_.end()) { int densityVal = cp->PointsWithinEps(); Carray densityArray; // Add all points of current density. while (cp->PointsWithinEps() == densityVal && cp != Points_.end()) densityArray.push_back( *(cp++) ); mprintf("Density value %i has %zu points.\n", densityVal, densityArray.size()); // Sort array by distance std::sort(densityArray.begin(), densityArray.end(), Cpoint::dist_sort()); // Take the average of the points weighted by their position. double wtDistAv = 0.0; double sumWts = 0.0; //std::vector<double> weights; //weights.reserve( densityArray.size() ); int maxPt = (int)densityArray.size() - 1; for (int ip = 0; ip != (int)densityArray.size(); ++ip) { double wt = exp( (double)(ip - maxPt) ); //mprintf("\t%10i %10u %10u %10g\n", densityVal, ip, maxPt, wt); wtDistAv += (densityArray[ip].Dist() * wt); sumWts += wt; //weights.push_back( wt ); } wtDistAv /= sumWts; // Calculate the weighted sample variance //double distSD = 0.0; //for (unsigned int ip = 0; ip != densityArray.size(); ++ip) { // double diff = densityArray[ip].Dist() - wtDistAv; // distSD += weights[ip] * (diff * diff); //} //distSD /= sumWts; weightedAverage.AddXY(densityVal, wtDistAv); //tempOut.Printf("%i %g %g %g\n", densityVal, wtDistAv, sqrt(distSD), sumWts); tempOut.Printf("%i %g %g\n", densityVal, wtDistAv, sumWts); /* // Find the median. double median, Q1, Q3; if (densityArray.size() == 1) { median = densityArray[0].Dist(); Q1 = median; Q3 = median; } else { unsigned int q3_beg; unsigned int med_idx = densityArray.size() / 2; // Always 0 <= Q1 < med_idx if ((densityArray.size() % 2) == 0) { median = (densityArray[med_idx].Dist() + densityArray[med_idx-1].Dist()) / 2.0; q3_beg = med_idx; } else { median = densityArray[med_idx].Dist(); q3_beg = med_idx + 1; } if (densityArray.size() == 2) { Q1 = densityArray[0].Dist(); Q3 = densityArray[1].Dist(); } else { // Find lower quartile unsigned int q1_idx = med_idx / 2; if ((med_idx % 2) == 0) Q1 = (densityArray[q1_idx].Dist() + densityArray[q1_idx-1].Dist()) / 2.0; else Q1 = densityArray[q1_idx].Dist(); // Find upper quartile unsigned int q3_size = densityArray.size() - q3_beg; unsigned int q3_idx = (q3_size / 2) + q3_beg; if ((q3_size %2) == 0) Q3 = (densityArray[q3_idx].Dist() + densityArray[q3_idx-1].Dist()) / 2.0; else Q3 = densityArray[q3_idx].Dist(); } } mprintf("\tMedian dist value is %g. Q1= %g Q3= %g\n", median, Q1, Q3); */ } tempOut.CloseFile(); // END CALCULATING WEIGHTED DISTANCE AVERAGE /* // TEST tempOut.OpenWrite("temp2.dat"); std::vector<double> Hist( Points_.back().PointsWithinEps()+1, 0.0 ); int gWidth = 3; double cval = 3.0; double two_c_squared = 2.0 * cval * cval; mprintf("DBG: cval= %g, Gaussian denominator is %g\n", cval, two_c_squared); for (int wtIdx = 0; wtIdx != (int)weightedAverage.Size(); wtIdx++) { int bval = weightedAverage.X(wtIdx); for (int xval = std::max(bval - gWidth, 0); xval != std::min(bval + gWidth + 1, (int)Hist.size()); xval++) { // a: height (weighted average) // b: center (density value) // c: width // x: density value in histogram //int xval = weightedAverage.X(idx); //double bval = weightedAverage.X(wtIdx); //double bval = (double)wtIdx; double diff = (double)(xval - bval); //Hist[xval] += (weightedAverage.Y(wtIdx) * exp( -( (diff * diff) / two_c_squared ) )); Hist[xval] = std::max(Hist[xval], weightedAverage.Y(wtIdx) * exp( -( (diff * diff) / two_c_squared ) )); } } for (unsigned int idx = 0; idx != Hist.size(); idx++) tempOut.Printf("%u %g\n", idx, Hist[idx]); tempOut.CloseFile(); // END TEST */ /* // TEST // Construct best-fit line segments tempOut.OpenWrite("temp2.dat"); double slope, intercept, correl; int segment_length = 3; DataSet_Mesh Segment; Segment.Allocate1D( segment_length ); for (int wtIdx = 0; wtIdx != (int)weightedAverage.Size(); wtIdx++) { Segment.Clear(); for (int idx = std::max(wtIdx - 1, 0); // TODO: use segment_length idx != std::min(wtIdx + 2, (int)weightedAverage.Size()); idx++) Segment.AddXY(weightedAverage.X(idx), weightedAverage.Y(idx)); Segment.LinearRegression(slope, intercept, correl, true); for (int idx = std::max(wtIdx - 1, 0); // TODO: use segment_length idx != std::min(wtIdx + 2, (int)weightedAverage.Size()); idx++) { double x = weightedAverage.X(idx); double y = slope * x + intercept; tempOut.Printf("%g %g %i\n", x, y, weightedAverage.X(wtIdx)); } } tempOut.CloseFile(); // END TEST */ // BEGIN WEIGHTED RUNNING AVG/SD OF DISTANCES // For each point, determine if it is greater than the average of the // weighted average distances of the previous, current, and next densities. int width = 2; int currentDensity = 0; int wtIdx = 0; double currentAvg = 0.0; double deltaSD = 0.0; double deltaAv = 0.0; int Ndelta = 0; CpptrajFile raOut; if (!rafile_.empty()) raOut.OpenWrite(rafile_); CpptrajFile raDelta; if (!radelta_.empty()) raDelta.OpenWrite(radelta_); std::vector<unsigned int> candidateIdxs; std::vector<double> candidateDeltas; cp = Points_.begin(); // Skip over points with zero density while (cp != Points_.end() && cp->PointsWithinEps() == 0) ++cp; while (weightedAverage.X(wtIdx) != cp->PointsWithinEps() && wtIdx < (int)Points_.size()) ++wtIdx; for (Carray::const_iterator point = cp; point != Points_.end(); ++point) { if (point->PointsWithinEps() != currentDensity) { //currentAvg = weightedAverage.Y(wtIdx); // New density value. Determine average. currentAvg = 0.0; // unsigned int Npt = 0; double currentWt = 0.0; for (int idx = std::max(wtIdx - width, 0); idx != std::min(wtIdx + width + 1, (int)weightedAverage.Size()); idx++) { //currentAvg += weightedAverage.Y(idx); //Npt++; double wt = weightedAverage.Y(idx); currentAvg += (weightedAverage.Y(idx) * wt); currentWt += wt; } //currentAvg /= (double)Npt; currentAvg /= currentWt; //smoothAv += currentAvg; //smoothSD += (currentAvg * currentAvg); //Nsmooth++; currentDensity = point->PointsWithinEps(); if (raOut.IsOpen()) raOut.Printf("%i %g %g\n", currentDensity, currentAvg, weightedAverage.Y(wtIdx)); wtIdx++; } double delta = (point->Dist() - currentAvg); if (delta > 0.0) { //delta *= log((double)currentDensity); if (raDelta.IsOpen()) raDelta.Printf("%8i %8.3f %8i %8.3f %8.3f\n", currentDensity, delta, point->Fnum()+1, point->Dist(), currentAvg); candidateIdxs.push_back( point - Points_.begin() ); candidateDeltas.push_back( delta ); deltaAv += delta; deltaSD += (delta * delta); Ndelta++; } } raOut.CloseFile(); deltaAv /= (double)Ndelta; deltaSD = (deltaSD / (double)Ndelta) - (deltaAv * deltaAv); if (deltaSD > 0.0) deltaSD = sqrt(deltaSD); else deltaSD = 0.0; if (raDelta.IsOpen()) raDelta.Printf("#DeltaAvg= %g DeltaSD= %g\n", deltaAv, deltaSD); raDelta.CloseFile(); int cnum = 0; for (unsigned int i = 0; i != candidateIdxs.size(); i++) { if (candidateDeltas[i] > (deltaSD)) { Points_[candidateIdxs[i]].SetCluster( cnum++ ); mprintf("\tPoint %u (frame %i, density %i) selected as candidate for cluster %i\n", candidateIdxs[i], Points_[candidateIdxs[i]].Fnum()+1, Points_[candidateIdxs[i]].PointsWithinEps(), cnum-1); } } // END WEIGHTED AVG/SD OF DISTANCES /* // Currently doing this by calculating the running average of density vs // distance, then choosing points with distance > twice the SD of the // running average. // NOTE: Store in a mesh data set for now in case we want to spline etc later. if (avg_factor_ < 1) avg_factor_ = 10; unsigned int window_size = Points_.size() / (unsigned int)avg_factor_; mprintf("\tRunning avg window size is %u\n", window_size); // FIXME: Handle case where window_size < frames DataSet_Mesh runavg; unsigned int ra_size = Points_.size() - window_size + 1; runavg.Allocate1D( ra_size ); double dwindow = (double)window_size; double sumx = 0.0; double sumy = 0.0; for (unsigned int i = 0; i < window_size; i++) { sumx += (double)Points_[i].PointsWithinEps(); sumy += Points_[i].Dist(); } runavg.AddXY( sumx / dwindow, sumy / dwindow ); for (unsigned int i = 1; i < ra_size; i++) { unsigned int nextwin = i + window_size - 1; unsigned int prevwin = i - 1; sumx = (double)Points_[nextwin].PointsWithinEps() - (double)Points_[prevwin].PointsWithinEps() + sumx; sumy = Points_[nextwin].Dist() - Points_[prevwin].Dist() + sumy; runavg.AddXY( sumx / dwindow, sumy / dwindow ); } // Write running average if (!rafile_.empty()) { CpptrajFile raOut; if (raOut.OpenWrite(rafile_)) mprinterr("Error: Could not open running avg file '%s' for write.\n", rafile_.c_str()); else { for (unsigned int i = 0; i != runavg.Size(); i++) raOut.Printf("%g %g\n", runavg.X(i), runavg.Y(i)); raOut.CloseFile(); } } double ra_sd; double ra_avg = runavg.Avg( ra_sd ); // Double stdev to use as cutoff for findning anomalously high peaks. ra_sd *= 2.0; mprintf("\tAvg of running avg set is %g, SD*2.0 (delta cutoff) is %g\n", ra_avg, ra_sd); // For each point in density vs distance plot, determine which running // average point is closest. If the difference between the point and the // running average point is > 2.0 the SD of the running average data, // consider it a 'peak'. CpptrajFile raDelta; if (!radelta_.empty()) raDelta.OpenWrite("radelta.dat"); if (raDelta.IsOpen()) raDelta.Printf("%-10s %10s %10s\n", "#Frame", "RnAvgPos", "Delta"); unsigned int ra_position = 0; // Position in the runavg DataSet unsigned int ra_end = runavg.Size() - 1; int cnum = 0; for (Carray::iterator point = Points_.begin(); point != Points_.end(); ++point) { if (ra_position != ra_end) { // Is the next running avgd point closer to this point? while (ra_position != ra_end) { double dens = (double)point->PointsWithinEps(); double diff0 = fabs( dens - runavg.X(ra_position ) ); double diff1 = fabs( dens - runavg.X(ra_position+1) ); if (diff1 < diff0) ++ra_position; // Next running avg position is closer for this point. else break; // This position is closer. } } double delta = point->Dist() - runavg.Y(ra_position); if (raDelta.IsOpen()) raDelta.Printf("%-10i %10u %10g", point->Fnum()+1, ra_position, delta); if (delta > ra_sd) { if (raDelta.IsOpen()) raDelta.Printf(" POTENTIAL CLUSTER %i", cnum); point->SetCluster(cnum++); } if (raDelta.IsOpen()) raDelta.Printf("\n"); } raDelta.CloseFile(); */ return cnum; }
// ----------------------------------------------------------------------------- int Cluster_DPeaks::Cluster_DiscreteDensity() { mprintf("\tStarting DPeaks clustering, discrete density calculation.\n"); Points_.clear(); // First determine which frames are being clustered. for (int frame = 0; frame < (int)FrameDistances_.Nframes(); ++frame) if (!FrameDistances_.IgnoringRow( frame )) Points_.push_back( Cpoint(frame) ); // Sanity check. if (Points_.size() < 2) { mprinterr("Error: Only 1 frame in initial clustering.\n"); return 1; } // For each point, determine how many others are within epsilon. Also // determine maximum distance between any two points. mprintf("\tDetermining local density of each point.\n"); ProgressBar cluster_progress( Points_.size() ); double maxDist = -1.0; for (Carray::iterator point0 = Points_.begin(); point0 != Points_.end(); ++point0) { cluster_progress.Update(point0 - Points_.begin()); int density = 0; for (Carray::const_iterator point1 = Points_.begin(); point1 != Points_.end(); ++point1) { if (point0 != point1) { double dist = FrameDistances_.GetFdist(point0->Fnum(), point1->Fnum()); maxDist = std::max(maxDist, dist); if ( dist < epsilon_ ) density++; } } point0->SetPointsWithinEps( density ); } mprintf("DBG: Max dist= %g\n", maxDist); // DEBUG: Frame/Density CpptrajFile fdout; fdout.OpenWrite("fd.dat"); for (Carray::const_iterator point = Points_.begin(); point != Points_.end(); ++point) fdout.Printf("%i %i\n", point->Fnum()+1, point->PointsWithinEps()); fdout.CloseFile(); // Sort by density here. Otherwise array indices will be invalid later. std::sort( Points_.begin(), Points_.end(), Cpoint::pointsWithinEps_sort() ); // For each point, find the closest point that has higher density. Since // array is now sorted by density the last point has the highest density. Points_.back().SetDist( maxDist ); mprintf("\tFinding closest neighbor point with higher density for each point.\n"); unsigned int lastidx = Points_.size() - 1; cluster_progress.SetupProgress( lastidx ); for (unsigned int idx0 = 0; idx0 != lastidx; idx0++) { cluster_progress.Update( idx0 ); double min_dist = maxDist; int nearestIdx = -1; // Index of nearest neighbor with higher density Cpoint& point0 = Points_[idx0]; //mprintf("\nDBG:\tSearching for nearest neighbor to idx %u with higher density than %i.\n", // idx0, point0.PointsWithinEps()); // Since array is sorted by density we can start at the next point. for (unsigned int idx1 = idx0+1; idx1 != Points_.size(); idx1++) { Cpoint const& point1 = Points_[idx1]; double dist1_2 = FrameDistances_.GetFdist(point0.Fnum(), point1.Fnum()); if (point1.PointsWithinEps() > point0.PointsWithinEps()) { if (dist1_2 < min_dist) { min_dist = dist1_2; nearestIdx = (int)idx1; //mprintf("DBG:\t\tNeighbor idx %i is closer (density %i, distance %g)\n", // nearestIdx, point1.PointsWithinEps(), min_dist); } } } point0.SetDist( min_dist ); //mprintf("DBG:\tClosest point to %u with higher density is %i (distance %g)\n", // idx0, nearestIdx, min_dist); point0.SetNearestIdx( nearestIdx ); } // Plot density vs distance for each point. if (!dvdfile_.empty()) { CpptrajFile output; if (output.OpenWrite(dvdfile_)) mprinterr("Error: Could not open density vs distance plot '%s' for write.\n", dvdfile_.c_str()); // TODO: Make fatal? else { output.Printf("%-10s %10s %s %10s %10s\n", "#Density", "Distance", "Frame", "Idx", "Neighbor"); for (Carray::const_iterator point = Points_.begin(); point != Points_.end(); ++point) output.Printf("%-10i %10g \"%i\" %10u %10i\n", point->PointsWithinEps(), point->Dist(), point->Fnum()+1, point-Points_.begin(), point->NearestIdx()); output.CloseFile(); } } return 0; }
// ----------------------------------------------------------------------------- int Cluster_DPeaks::Cluster_GaussianKernel() { mprintf("\tStarting DPeaks clustering. Using Gaussian kernel to calculate density.\n"); // First determine which frames are being clustered. Points_.clear(); int oidx = 0; for (int frame = 0; frame < (int)FrameDistances_.Nframes(); ++frame) if (!FrameDistances_.IgnoringRow( frame )) Points_.push_back( Cpoint(frame, oidx++) ); // Sanity check. if (Points_.size() < 2) { mprinterr("Error: Only 1 frame in initial clustering.\n"); return 1; } // Sort distances std::vector<float> Distances; for (ClusterMatrix::const_iterator mat = FrameDistances_.begin(); mat != FrameDistances_.end(); ++mat) Distances.push_back( *mat ); std::sort( Distances.begin(), Distances.end() ); unsigned int idx = (unsigned int)((double)Distances.size() * 0.02); double bandwidth = (double)Distances[idx]; mprintf("idx= %u, bandwidth= %g\n", idx, bandwidth); // Density via Gaussian kernel double maxDist = -1.0; for (unsigned int i = 0; i != Points_.size(); i++) { for (unsigned int j = i+1; j != Points_.size(); j++) { double dist = FrameDistances_.GetFdist(Points_[i].Fnum(), Points_[j].Fnum()); maxDist = std::max( maxDist, dist ); dist /= bandwidth; double gk = exp(-(dist *dist)); Points_[i].AddDensity( gk ); Points_[j].AddDensity( gk ); } } mprintf("Max dist= %g\n", maxDist); CpptrajFile rhoOut; rhoOut.OpenWrite("rho.dat"); for (unsigned int i = 0; i != Points_.size(); i++) rhoOut.Printf("%u %g\n", i+1, Points_[i].Density()); rhoOut.CloseFile(); // Sort by density, descending std::stable_sort( Points_.begin(), Points_.end(), Cpoint::density_sort_descend() ); CpptrajFile ordrhoOut; ordrhoOut.OpenWrite("ordrho.dat"); for (unsigned int i = 0; i != Points_.size(); i++) ordrhoOut.Printf("%u %g %i %i\n", i+1, Points_[i].Density(), Points_[i].Fnum()+1, Points_[i].Oidx()+1); ordrhoOut.CloseFile(); // Determine minimum distances int first_idx = Points_[0].Oidx(); Points_[first_idx].SetDist( -1.0 ); Points_[first_idx].SetNearestIdx(-1); for (unsigned int ii = 1; ii != Points_.size(); ii++) { int ord_i = Points_[ii].Oidx(); Points_[ord_i].SetDist( maxDist ); for (unsigned int jj = 0; jj != ii; jj++) { int ord_j = Points_[jj].Oidx(); double dist = FrameDistances_.GetFdist(Points_[ord_i].Fnum(), Points_[ord_j].Fnum()); if (dist < Points_[ord_i].Dist()) { Points_[ord_i].SetDist( dist ); Points_[ord_j].SetNearestIdx( ord_j ); } } } if (!dvdfile_.empty()) { CpptrajFile output; if (output.OpenWrite(dvdfile_)) return 1; for (Carray::const_iterator point = Points_.begin(); point != Points_.end(); ++point) output.Printf("%g %g %i\n", point->Density(), point->Dist(), point->NearestIdx()+1); output.CloseFile(); } return 0; }
/** Given the structure of a molecule and its normal mode vibrational * frequencies this routine uses standard statistical mechanical * formulas for an ideal gas (in the canonical ensemble, see, * for example, d. a. mcquarrie, "statistical thermodynamics", * harper & row, new york, 1973, chapters 5, 6, and 8) to compute * the entropy, heat capacity, and internal energy. * The si system of units is used internally. Conversion to units * more familiar to most chemists is made for output. * * \param outfile output file, should already be open. * \param natoms Number of atoms * \param nvecs Number of eigenvectors (already converted to frequencies) * \param crd coordinates in Angstroms * \param amass atomic weights, in amu. * \param freq vibrational frequencies, in cm**-1 and in ascending order * \param temp temperature * \param patm pressure, in atmospheres */ void thermo( CpptrajFile& outfile, int natoms, int nvecs, int ilevel, const double* crd, const double* amass, const double* freq, double temp, double patm) { // pmom principal moments of inertia, in amu-bohr**2 and in ascending order. double pmom[3], rtemp, rtemp1, rtemp2, rtemp3; // ----- Constants ------------------- const double thresh = 900.0; // vibrational frequency threshold const double tokg = 1.660531e-27; // kilograms per amu. const double boltz = 1.380622e-23; // boltzman constant, in joules per kelvin. const double planck = 6.626196e-34; // planck constant, in joule-seconds. const double avog = 6.022169e+23; // avogadro constant, in mol**(-1). const double jpcal = 4.18674e+00; // joules per calorie. const double tomet = 1.0e-10; // metres per Angstrom. const double hartre = 4.35981e-18; // joules per hartree. const double pstd = 1.01325e+05; // Standard pressure in pascals // ----------------------------------- // compute the gas constant, pi, pi**2, and e. // compute the conversion factors cal per joule and kcal per joule. const double gas = avog * boltz; // pi = four * datan(one) const double pipi = PI * PI; const double e = exp(1.0); const double tocal = 1.0 / jpcal; const double tokcal = tocal / 1000.0; if (!outfile.IsOpen()) { mprinterr("Internal Error: thermo: output file is not open.\n"); return; } // print the temperature and pressure. outfile.Printf("\n *******************\n"); outfile.Printf( " - Thermochemistry -\n"); outfile.Printf( " *******************\n\n"); outfile.Printf("\n temperature %9.3f kelvin\n pressure %9.5f atm\n",temp,patm); double pressure = pstd * patm; double rt = gas * temp; // compute and print the molecular mass in amu, then convert to // kilograms. double weight = 0.0; for (int iat = 0; iat < natoms; ++iat) weight += amass[iat]; outfile.Printf(" molecular mass (principal isotopes) %11.5f amu\n", weight); weight *= tokg; //trap non-unit multiplicities. //if (multip != 1) { // outfile.Printf("\n Warning-- assumptions made about the electronic partition function\n"); // outfile.Printf( " are not valid for multiplets!\n\n"); //} // compute contributions due to translation: // etran-- internal energy // ctran-- constant v heat capacity // stran-- entropy double dum1 = boltz * temp; double dum2 = pow(TWOPI, 1.5); double arg = pow(dum1, 1.5) / planck; arg = (arg / pressure) * (dum1 / planck); arg = arg * dum2 * (weight / planck); arg = arg * sqrt(weight) * exp(2.5); double stran = gas * log(arg); double etran = 1.5 * rt; double ctran = 1.5 * gas; // Compute contributions due to electronic motion: // It is assumed that the first electronic excitation energy // is much greater than kt and that the ground state has a // degeneracy of one. Under these conditions the electronic // partition function can be considered to be unity. The // ground electronic state is taken to be the zero of // electronic energy. // for monatomics print and return. if (natoms <= 1){ outfile.Printf("\n internal energy: %10.3f joule/mol %10.3f kcal/mol\n", etran, etran * tokcal); outfile.Printf( " entropy: %10.3f joule/k-mol %10.3f cal/k-mol\n", stran, stran * tocal); outfile.Printf( " heat capacity cv: %10.3f joule/k-mol %10.3f cal/k-mol\n", ctran, ctran * tocal); return; } // Allocate workspace memory // vtemp vibrational temperatures, in kelvin. // evibn contribution to e from the vibration n. // cvibn contribution to cv from the vibration n. // svibn contribution to s from the vibration n. double* WorkSpace = new double[ 4 * nvecs ]; double* vtemp = WorkSpace; double* evibn = WorkSpace + nvecs; double* cvibn = WorkSpace + nvecs*2; double* svibn = WorkSpace + nvecs*3; // compute contributions due to rotation. // Compute the principal moments of inertia, get the rotational // symmetry number, see if the molecule is linear, and compute // the rotational temperatures. Note the imbedded conversion // of the moments to SI units. MomentOfInertia( natoms, crd, amass, pmom ); outfile.Printf("\n principal moments of inertia (nuclei only) in amu-A**2:\n"); outfile.Printf( " %12.2f%12.2f%12.2f\n", pmom[0], pmom[1], pmom[2]); bool linear = false; // Symmetry number: only for linear molecules. for others symmetry number is unity double sn = 1.0; if (natoms <= 2) { linear = true; if (amass[0] == amass[1]) sn = 2.0; } outfile.Printf("\n rotational symmetry number %3.0f\n", sn); double con = planck / (boltz*8.0*pipi); con = (con / tokg) * (planck / (tomet*tomet)); if (linear) { rtemp = con / pmom[2]; if (rtemp < 0.2) { outfile.Printf("\n Warning-- assumption of classical behavior for rotation\n"); outfile.Printf( " may cause significant error\n"); } outfile.Printf("\n rotational temperature (kelvin) %12.5f\n", rtemp); } else { rtemp1 = con / pmom[0]; rtemp2 = con / pmom[1]; rtemp3 = con / pmom[2]; if (rtemp1 < 0.2) { outfile.Printf("\n Warning-- assumption of classical behavior for rotation\n"); outfile.Printf( " may cause significant error\n"); } outfile.Printf("\n rotational temperatures (kelvin) %12.5f%12.5f%12.5f\n", rtemp1, rtemp2, rtemp3); } // erot-- rotational contribution to internal energy. // crot-- rotational contribution to cv. // srot-- rotational contribution to entropy. double erot, crot, srot; if (linear) { erot = rt; crot = gas; arg = (temp/rtemp) * (e/sn); srot = gas * log(arg); } else { erot = 1.5 * rt; crot = 1.5 * gas; arg = sqrt(PI*e*e*e) / sn; double dum = (temp/rtemp1) * (temp/rtemp2) * (temp/rtemp3); arg = arg * sqrt(dum); srot = gas * log(arg); } // compute contributions due to vibration. // compute vibrational temperatures and zero point vibrational // energy. only real frequencies are included in the analysis. // ndof = 3*natoms - 6 - nimag // if (nimag .ne. 0) write(iout,1210) nimag // if (linear) ndof = ndof + 1 int ndof = nvecs; // (---iff is the first frequency to include in thermo:) int iff; if (ilevel != 0) iff = 0; else if (linear) iff = 5; else iff = 6; con = planck / boltz; double ezpe = 0.0; for (int i = 0; i < ndof; ++i) { vtemp[i] = freq[i+iff] * con * 3.0e10; ezpe += freq[i+iff] * 3.0e10; } ezpe = 0.5 * planck * ezpe; outfile.Printf("\n zero point vibrational energy %12.1f (joules/mol) \n",ezpe * avog); outfile.Printf( " %12.5f (kcal/mol)\n",ezpe * tokcal * avog); outfile.Printf( " %12.7f (hartree/particle)\n", ezpe / hartre); // compute the number of vibrations for which more than 5% of an // assembly of molecules would exist in vibrational excited states. // special printing for these modes is done to allow the user to // easily take internal rotations into account. the criterion // corresponds roughly to a low frequency of 1.9(10**13) hz, or // 625 cm**(-1), or a vibrational temperature of 900 k. int lofreq = 0; for (int i = 0; i < ndof; ++i) if (vtemp[i] < thresh) ++lofreq; if (lofreq != 0) { outfile.Printf("\n Warning-- %3i vibrations have low frequencies and may represent hindered \n", lofreq); outfile.Printf( " internal rotations. The contributions printed below assume that these \n"); outfile.Printf( " really are vibrations.\n"); } // compute: // evib-- the vibrational component of the internal energy. // cvib-- the vibrational component of the heat capacity. // svib-- the vibrational component of the entropy. double evib = 0.0; double cvib = 0.0; double svib = 0.0; double scont; for (int i = 0; i < ndof; ++i) { // compute some common factors. double tovt = vtemp[i] / temp; double etovt = exp(tovt); double em1 = etovt - 1.0; // compute contributions due to the i'th vibration. double econt = tovt * (0.5 + 1.0/em1); double ccont = etovt * pow(tovt/em1,2.0); double argd = 1.0 - 1.0/etovt; if (argd > 1.0e-7) scont = tovt/em1 - log(argd); else { scont = 0.0; outfile.Printf(" warning: setting vibrational entropy to zero for mode %i with vtemp = %f\n", i+1, vtemp[i]); } // if (lofreq .ge. i) then evibn[i] = econt * rt; cvibn[i] = ccont * gas; svibn[i] = scont * gas; // end if evib += econt; cvib += ccont; svib += scont; } evib *= rt; cvib *= gas; svib *= gas; // the units are now: // e-- joules/mol // c-- joules/mol-kelvin // s-- joules/mol-kelvin double etot = etran + erot + evib; double ctot = ctran + crot + cvib; double stot = stran + srot + svib; // print the sum of the hartree-fock energy and the thermal energy. // call tread(501,gen,47,1,47,1,0) // esum = gen(32) + etot/avog/hartre // write(iout,1230) esum // convert to the following and print // e-- kcal/mol // c-- cal/mol-kelvin // s-- cal/mol-kelvin etran = etran * tokcal; ctran = ctran * tocal; stran = stran * tocal; erot = erot * tokcal; crot = crot * tocal; srot = srot * tocal; evib = evib * tokcal; cvib = cvib * tocal; svib = svib * tocal; etot = etran + erot + evib; ctot = ctran + crot + cvib; stot = stran + srot + svib; for (int i = 0; i < ndof; ++i) { evibn[i] *= tokcal; cvibn[i] *= tocal; svibn[i] *= tocal; } outfile.Printf("\n\n freq. E Cv S\n"); outfile.Printf( " cm**-1 kcal/mol cal/mol-kelvin cal/mol-kelvin\n"); outfile.Printf( "--------------------------------------------------------------------------------\n"); outfile.Printf( " Total %11.3f %11.3f %11.3f\n",etot,ctot,stot); outfile.Printf( " translational %11.3f %11.3f %11.3f\n",etran,ctran,stran); outfile.Printf( " rotational %11.3f %11.3f %11.3f\n",erot,crot,srot); outfile.Printf( " vibrational %11.3f %11.3f %11.3f\n",evib,cvib,svib); for (int i = 0; i < iff; ++i) outfile.Printf(" %5i%10.3f\n", i+1, freq[i]); for (int i = 0; i < ndof; ++i) { outfile.Printf(" %5i%10.3f %11.3f %11.3f %11.3f\n",i+iff+1, freq[i+iff], evibn[i], cvibn[i], svibn[i]); } delete[] WorkSpace; }
// Action_AtomMap::Init() Action::RetType Action_AtomMap::Init(ArgList& actionArgs, ActionInit& init, int debugIn) { DataFile* rmsout = 0; int refatom,targetatom; debug_ = debugIn; RefMap_.SetDebug(debug_); TgtMap_.SetDebug(debug_); // Get Args CpptrajFile* outputfile = init.DFL().AddCpptrajFile(actionArgs.GetStringKey("mapout"), "Atom Map"); maponly_ = actionArgs.hasKey("maponly"); rmsfit_ = actionArgs.hasKey("rmsfit"); if (rmsfit_) rmsout = init.DFL().AddDataFile( actionArgs.GetStringKey("rmsout"), actionArgs ); std::string targetName = actionArgs.GetStringNext(); std::string refName = actionArgs.GetStringNext(); if (targetName.empty()) { mprinterr("Error: No target specified.\n"); return Action::ERR; } if (refName.empty()) { mprinterr("Error: No reference specified.\n"); return Action::ERR; } // Get Reference RefFrame_ = (DataSet_Coords_REF*)init.DSL().FindSetOfType( refName, DataSet::REF_FRAME ); if (RefFrame_ == 0) { mprinterr("Error: Could not get reference frame %s\n",refName.c_str()); return Action::ERR; } // Get Target TgtFrame_ = (DataSet_Coords_REF*)init.DSL().FindSetOfType( targetName, DataSet::REF_FRAME ); if (TgtFrame_ == 0) { mprinterr("Error: Could not get target frame %s\n",targetName.c_str()); return Action::ERR; } mprintf(" ATOMMAP: Atoms in trajectories associated with parm %s will be\n", TgtFrame_->Top().c_str()); mprintf(" mapped according to parm %s.\n",RefFrame_->Top().c_str()); if (outputfile != 0) mprintf(" Map will be written to %s\n",outputfile->Filename().full()); if (maponly_) mprintf(" maponly: Map will only be written, not used in trajectory read.\n"); if (!maponly_ && rmsfit_) { mprintf(" rmsfit: Will rms fit mapped atoms in tgt to reference.\n"); if (rmsout != 0) { rmsdata_ = init.DSL().AddSet(DataSet::DOUBLE, actionArgs.GetStringNext(), "RMSD"); if (rmsdata_==0) return Action::ERR; rmsout->AddDataSet(rmsdata_); } } // For each map, set up (get element for each atom, initialize map mem), // determine what atoms are bonded to each other via simple distance // cutoffs, the give each atom an ID based on what atoms are bonded to // it, noting which IDs are unique for that map. if (RefMap_.Setup(RefFrame_->Top())!=0) return Action::ERR; //RefMap_.WriteMol2((char*)"RefMap.mol2\0"); // DEBUG RefMap_.DetermineAtomIDs(); if (TgtMap_.Setup(TgtFrame_->Top())!=0) return Action::ERR; //TgtMap_.WriteMol2((char*)"TgtMap.mol2\0"); // DEBUG TgtMap_.DetermineAtomIDs(); // Check if number of atoms in each map is equal if (RefMap_.Natom() != TgtMap_.Natom()) { mprintf("Warning: # atoms in reference (%i) not equal\n", RefMap_.Natom()); mprintf("Warning:\tto # atoms in target (%i).\n",TgtMap_.Natom()); } // Set up RMS frames to be able to hold max # of possible atoms rmsRefFrame_.SetupFrame(RefMap_.Natom()); rmsTgtFrame_.SetupFrame(RefMap_.Natom()); // Allocate memory for atom map // AMap_[reference]=target AMap_.resize( RefMap_.Natom(), -1); // Map unique atoms int numMappedAtoms = MapUniqueAtoms(RefMap_, TgtMap_); if (debug_>0) mprintf("* MapUniqueAtoms: %i atoms mapped.\n",numMappedAtoms); // If no unique atoms mapped system is highly symmetric and needs to be // iteratively mapped. Otherwise just map remaining atoms. if (numMappedAtoms==0) { if (MapWithNoUniqueAtoms(RefMap_,TgtMap_)) return Action::ERR; } else { if (MapAtoms(RefMap_,TgtMap_)) return Action::ERR; } // Print atom map and count # mapped atoms numMappedAtoms = 0; outputfile->Printf("%-6s %4s %6s %4s\n","#TgtAt","Tgt","RefAt","Ref"); for (refatom=0; refatom < RefMap_.Natom(); refatom++) { targetatom = AMap_[refatom]; if (targetatom < 0) outputfile->Printf("%6s %4s %6i %4s\n","---","---",refatom+1,RefMap_[refatom].c_str()); else outputfile->Printf("%6i %4s %6i %4s\n",targetatom+1,TgtMap_[targetatom].c_str(), refatom+1, RefMap_[refatom].c_str()); if (targetatom>=0) { //mprintf("* TargetAtom %6i(%4s) maps to RefAtom %6i(%4s)\n", // targetatom,TgtMap_.P->names[targetatom], // refatom,RefMap_.P->names[refatom]); ++numMappedAtoms; } //else { // mprintf("* Could not map any TargetAtom to RefAtom %6i(%4s)\n", // refatom,RefMap_.P->names[refatom]); //} } mprintf(" %i total atoms were mapped.\n",numMappedAtoms); if (maponly_) return Action::OK; // If rmsfit is specified, an rms fit of target to reference will be // performed using all atoms that were successfully mapped from // target to reference. if (rmsfit_) { // Set up a reference frame containing only mapped reference atoms rmsRefFrame_.StripUnmappedAtoms(RefFrame_->RefFrame(), AMap_); mprintf(" rmsfit: Will rms fit %i atoms from target to reference.\n",numMappedAtoms); return Action::OK; } // Check if not all atoms could be mapped if (numMappedAtoms != RefMap_.Natom()) { // If the number of mapped atoms is less than the number of reference // atoms but equal to the number of target atoms, can modify the reference // frame to only include mapped atoms if (numMappedAtoms<RefMap_.Natom() && numMappedAtoms==TgtMap_.Natom()) { // Create mask that includes only reference atoms that could be mapped AtomMask M1; for (refatom = 0; refatom < RefMap_.Natom(); refatom++) { if (AMap_[refatom] != -1) M1.AddAtom(refatom); } // Strip reference parm mprintf(" Modifying reference '%s' topology and frame to match mapped atoms.\n", RefFrame_->legend()); if (RefFrame_->StripRef( M1 )) return Action::ERR; // Since AMap[ ref ] = tgt but ref is now missing any stripped atoms, // the indices of AMap must be shifted to match int refIndex = 0; // The new index for (refatom = 0; refatom < RefMap_.Natom(); refatom++) { targetatom = AMap_[refatom]; if (targetatom<0) continue; else AMap_[refIndex++]=targetatom; } } else { mprintf("Warning: AtomMap: Not all atoms were mapped. Frames will not be modified.\n"); maponly_=true; } } if (!maponly_) { // Set up new Frame newFrame_ = new Frame(); newFrame_->SetupFrameM( TgtFrame_->Top().Atoms() ); // Set up new Parm newParm_ = TgtFrame_->Top().ModifyByMap(AMap_); } return Action::OK; }
int Parm_CharmmPsf::WriteParm(FileName const& fname, Topology const& parm) { // TODO: CMAP etc info CpptrajFile outfile; if (outfile.OpenWrite(fname)) return 1; // Write PSF outfile.Printf("PSF\n\n"); // Write title std::string titleOut = parm.ParmName(); titleOut.resize(78); outfile.Printf("%8i !NTITLE\n* %-78s\n\n", 1, titleOut.c_str()); // Write NATOM section outfile.Printf("%8i !NATOM\n", parm.Natom()); unsigned int idx = 1; // Make fake segment ids for now. char segid[2]; segid[0] = 'A'; segid[1] = '\0'; mprintf("Warning: Assigning single letter segment IDs.\n"); int currentMol = 0; bool inSolvent = false; for (Topology::atom_iterator atom = parm.begin(); atom != parm.end(); ++atom, ++idx) { int resnum = atom->ResNum(); if (atom->MolNum() != currentMol) { if (!inSolvent) { inSolvent = parm.Mol(atom->MolNum()).IsSolvent(); currentMol = atom->MolNum(); segid[0]++; } else inSolvent = parm.Mol(atom->MolNum()).IsSolvent(); } // TODO: Print type name for xplor-like PSF int typeindex = atom->TypeIndex() + 1; // If type begins with digit, assume charmm numbers were read as // type. Currently Amber types all begin with letters. if (isdigit(atom->Type()[0])) typeindex = convertToInteger( *(atom->Type()) ); // ATOM# SEGID RES# RES ATNAME ATTYPE CHRG MASS (REST OF COLUMNS ARE LIKELY FOR CMAP AND CHEQ) outfile.Printf("%8i %-4s %-4i %-4s %-4s %4i %14.6G %9g %10i\n", idx, segid, parm.Res(resnum).OriginalResNum(), parm.Res(resnum).c_str(), atom->c_str(), typeindex, atom->Charge(), atom->Mass(), 0); } outfile.Printf("\n"); // Write NBOND section outfile.Printf("%8u !NBOND: bonds\n", parm.Bonds().size() + parm.BondsH().size()); idx = 1; for (BondArray::const_iterator bond = parm.BondsH().begin(); bond != parm.BondsH().end(); ++bond, ++idx) { outfile.Printf("%8i%8i", bond->A1()+1, bond->A2()+1); if ((idx % 4)==0) outfile.Printf("\n"); } for (BondArray::const_iterator bond = parm.Bonds().begin(); bond != parm.Bonds().end(); ++bond, ++idx) { outfile.Printf("%8i%8i", bond->A1()+1, bond->A2()+1); if ((idx % 4)==0) outfile.Printf("\n"); } if ((idx % 4)!=0) outfile.Printf("\n"); outfile.Printf("\n"); // Write NTHETA section outfile.Printf("%8u !NTHETA: angles\n", parm.Angles().size() + parm.AnglesH().size()); idx = 1; for (AngleArray::const_iterator ang = parm.AnglesH().begin(); ang != parm.AnglesH().end(); ++ang, ++idx) { outfile.Printf("%8i%8i%8i", ang->A1()+1, ang->A2()+1, ang->A3()+1); if ((idx % 3)==0) outfile.Printf("\n"); } for (AngleArray::const_iterator ang = parm.Angles().begin(); ang != parm.Angles().end(); ++ang, ++idx) { outfile.Printf("%8i%8i%8i", ang->A1()+1, ang->A2()+1, ang->A3()+1); if ((idx % 3)==0) outfile.Printf("\n"); } if ((idx % 3)==0) outfile.Printf("\n"); outfile.Printf("\n"); // Write out NPHI section outfile.Printf("%8u !NPHI: dihedrals\n", parm.Dihedrals().size() + parm.DihedralsH().size()); idx = 1; for (DihedralArray::const_iterator dih = parm.DihedralsH().begin(); dih != parm.DihedralsH().end(); ++dih, ++idx) { outfile.Printf("%8i%8i%8i%8i", dih->A1()+1, dih->A2()+1, dih->A3()+1, dih->A4()+1); if ((idx % 2)==0) outfile.Printf("\n"); } for (DihedralArray::const_iterator dih = parm.Dihedrals().begin(); dih != parm.Dihedrals().end(); ++dih, ++idx) { outfile.Printf("%8i%8i%8i%8i", dih->A1()+1, dih->A2()+1, dih->A3()+1, dih->A4()+1); if ((idx % 2)==0) outfile.Printf("\n"); } if ((idx % 2)==0) outfile.Printf("\n"); outfile.Printf("\n"); outfile.CloseFile(); return 0; }
// DataIO_Std::WriteDataNormal() int DataIO_Std::WriteDataNormal(CpptrajFile& file, DataSetList const& Sets) { // Assume all 1D data sets. if (Sets.empty() || CheckAllDims(Sets, 1)) return 1; // For this output to work the X-dimension of all sets needs to match. // The most important things for output are min and step so just check that. // Use X dimension of set 0 for all set dimensions. CheckXDimension( Sets ); // TODO: Check for empty dim. DataSet* Xdata = Sets[0]; Dimension const& Xdim = static_cast<Dimension const&>( Xdata->Dim(0) ); int xcol_width = Xdim.Label().size() + 1; // Only used if hasXcolumn_ if (xcol_width < 8) xcol_width = 8; int xcol_precision = 3; // Determine size of largest DataSet. size_t maxFrames = DetermineMax( Sets ); // Set up X column. TextFormat x_col_format(XcolFmt()); if (hasXcolumn_) { if (XcolPrecSet()) { xcol_width = XcolWidth(); x_col_format = TextFormat(XcolFmt(), XcolWidth(), XcolPrec()); } else { // Create format string for X column based on dimension in first data set. // Adjust X col precision as follows: if the step is set and has a // fractional component set the X col width/precision to either the data // width/precision or the current width/precision, whichever is larger. If // the set is XYMESH but step has not been set (so we do not know spacing // between X values) use default precision. Otherwise the step has no // fractional component so make the precision zero. double step_i; double step_f = modf( Xdim.Step(), &step_i ); double min_f = modf( Xdim.Min(), &step_i ); if (Xdim.Step() > 0.0 && (step_f > 0.0 || min_f > 0.0)) { xcol_precision = std::max(xcol_precision, Xdata->Format().Precision()); xcol_width = std::max(xcol_width, Xdata->Format().Width()); } else if (Xdata->Type() != DataSet::XYMESH) xcol_precision = 0; x_col_format.SetCoordFormat( maxFrames, Xdim.Min(), Xdim.Step(), xcol_width, xcol_precision ); } } else { // If not writing an X-column, no leading space for the first dataset. Xdata->SetupFormat().SetFormatAlign( TextFormat::RIGHT ); } // Write header to buffer std::vector<int> Original_Column_Widths; if (writeHeader_) { // If x-column present, write x-label if (hasXcolumn_) WriteNameToBuffer( file, Xdim.Label(), xcol_width, true ); // To prevent truncation of DataSet legends, adjust the width of each // DataSet if necessary. for (DataSetList::const_iterator ds = Sets.begin(); ds != Sets.end(); ++ds) { // Record original column widths in case they are changed. Original_Column_Widths.push_back( (*ds)->Format().Width() ); int colLabelSize; if (ds == Sets.begin() && !hasXcolumn_) colLabelSize = (int)(*ds)->Meta().Legend().size() + 1; else colLabelSize = (int)(*ds)->Meta().Legend().size(); //mprintf("DEBUG: Set '%s', fmt width= %i, colWidth= %i, colLabelSize= %i\n", // (*ds)->legend(), (*ds)->Format().Width(), (*ds)->Format().ColumnWidth(), // colLabelSize); if (colLabelSize >= (*ds)->Format().ColumnWidth()) (*ds)->SetupFormat().SetFormatWidth( colLabelSize ); } // Write dataset names to header, left-aligning first set if no X-column DataSetList::const_iterator set = Sets.begin(); if (!hasXcolumn_) WriteNameToBuffer( file, (*set)->Meta().Legend(), (*set)->Format().ColumnWidth(), true ); else WriteNameToBuffer( file, (*set)->Meta().Legend(), (*set)->Format().ColumnWidth(), false ); ++set; for (; set != Sets.end(); ++set) WriteNameToBuffer( file, (*set)->Meta().Legend(), (*set)->Format().ColumnWidth(), false ); file.Printf("\n"); } // Write Data DataSet::SizeArray positions(1); for (positions[0] = 0; positions[0] < maxFrames; positions[0]++) { // Output Frame for each set if (hasXcolumn_) file.Printf( x_col_format.fmt(), Xdata->Coord(0, positions[0]) ); for (DataSetList::const_iterator set = Sets.begin(); set != Sets.end(); ++set) (*set)->WriteBuffer(file, positions); file.Printf("\n"); } // Restore original column widths if necessary if (!Original_Column_Widths.empty()) for (unsigned int i = 0; i != Original_Column_Widths.size(); i++) Sets[i]->SetupFormat().SetFormatWidth( Original_Column_Widths[i] ); return 0; }
void Cluster_ReadInfo::ClusterResults(CpptrajFile& outfile) const { outfile.Printf("#Algorithm: Read from file '%s'\n", filename_.c_str()); if (!algorithm_.empty()) outfile.Printf("#Original algorithm: %s", algorithm_.c_str()); }
// Cluster_DBSCAN::ComputeKdistMap() void Cluster_DBSCAN::ComputeKdistMap( Range const& Kvals, std::vector<int> const& FramesToCluster ) const { int pt1_idx, pt2_idx, d_idx, point; mprintf("\tCalculating Kdist map for %s\n", Kvals.RangeArg()); double* kdist_array; // Store distance of pt1 to every other point. int nframes = (int)FramesToCluster.size(); // Ensure all Kdist points are within proper range Range::const_iterator kval; for (kval = Kvals.begin(); kval != Kvals.end(); ++kval) if (*kval < 1 || *kval >= nframes) { mprinterr("Error: Kdist value %i is out of range (1 <= Kdist < %i)\n", *kval, nframes); return; } int nvals = (int)Kvals.Size(); double** KMAP; // KMAP[i] has the ith nearest point for each point. KMAP = new double*[ nvals ]; for (int i = 0; i != nvals; i++) KMAP[i] = new double[ nframes ]; ParallelProgress progress( nframes ); # ifdef _OPENMP # pragma omp parallel private(pt1_idx, pt2_idx, d_idx, kval, point, kdist_array) firstprivate(progress) { progress.SetThread( omp_get_thread_num() ); #endif kdist_array = new double[ nframes ]; # ifdef _OPENMP # pragma omp for # endif for (pt1_idx = 0; pt1_idx < nframes; pt1_idx++) // X { progress.Update( pt1_idx ); point = FramesToCluster[pt1_idx]; d_idx = 0; // Store distances from pt1 to pt2 for (pt2_idx = 0; pt2_idx != nframes; pt2_idx++) kdist_array[d_idx++] = FrameDistances_.GetFdist(point, FramesToCluster[pt2_idx]); // Sort distances; will be smallest to largest std::sort( kdist_array, kdist_array + nframes ); // Save the distance of specified nearest neighbors to this point. d_idx = 0; for (kval = Kvals.begin(); kval != Kvals.end(); ++kval) // Y KMAP[d_idx++][pt1_idx] = kdist_array[ *kval ]; } delete[] kdist_array; # ifdef _OPENMP } // END omp parallel # endif progress.Finish(); // Sort all of the individual kdist plots, smallest to largest. for (int i = 0; i != nvals; i++) std::sort(KMAP[i], KMAP[i] + nframes); // Save in matrix, largest to smallest. DataSet_MatrixDbl kmatrix; kmatrix.Allocate2D( FramesToCluster.size(), Kvals.Size() ); for (int y = 0; y != nvals; y++) { for (int x = nframes - 1; x != -1; x--) kmatrix.AddElement( KMAP[y][x] ); delete[] KMAP[y]; } delete[] KMAP; // Write matrix to file DataFile outfile; ArgList outargs("usemap"); outfile.SetupDatafile(k_prefix_ + "Kmatrix.gnu", outargs, debug_); outfile.AddDataSet( (DataSet*)&kmatrix ); outfile.WriteDataOut(); // Write out the largest and smallest values for each K. // This means for each value of K the point with the furthest Kth-nearest // neighbor etc. CpptrajFile maxfile; if (maxfile.OpenWrite(k_prefix_ + "Kmatrix.max.dat")) return; maxfile.Printf("%-12s %12s %12s\n", "#Kval", "MaxD", "MinD"); d_idx = 0; for (kval = Kvals.begin(); kval != Kvals.end(); ++kval, d_idx++) maxfile.Printf("%12i %12g %12g\n", *kval, kmatrix.GetElement(0, d_idx), kmatrix.GetElement(nframes-1, d_idx)); maxfile.CloseFile(); }
// DataIO_Std::WriteSet3D() int DataIO_Std::WriteSet3D( DataSet const& setIn, CpptrajFile& file ) { if (setIn.Ndim() != 3) { mprinterr("Internal Error: DataSet %s in DataFile %s has %zu dimensions, expected 3.\n", setIn.legend(), file.Filename().full(), setIn.Ndim()); return 1; } DataSet_3D const& set = static_cast<DataSet_3D const&>( setIn ); Dimension const& Xdim = static_cast<Dimension const&>(set.Dim(0)); Dimension const& Ydim = static_cast<Dimension const&>(set.Dim(1)); Dimension const& Zdim = static_cast<Dimension const&>(set.Dim(2)); //if (Xdim.Step() == 1.0) xcol_precision = 0; if (sparse_) mprintf("\tOnly writing voxels with value > %g\n", cut_); // Print X Y Z Values // x y z val(x,y,z) DataSet::SizeArray pos(3); if (writeHeader_) { file.Printf("#counts %zu %zu %zu\n", set.NX(), set.NY(), set.NZ()); file.Printf("#origin %12.7f %12.7f %12.7f\n", set.Bin().GridOrigin()[0], set.Bin().GridOrigin()[1], set.Bin().GridOrigin()[2]); if (set.Bin().IsOrthoGrid()) { GridBin_Ortho const& b = static_cast<GridBin_Ortho const&>( set.Bin() ); file.Printf("#delta %12.7f %12.7f %12.7f\n", b.DX(), b.DY(), b.DZ()); } else { GridBin_Nonortho const& b = static_cast<GridBin_Nonortho const&>( set.Bin() ); file.Printf("#delta %12.7f %12.7f %12.7f %12.7f %12.7f %12.7f %12.7f %12.7f %12.7f\n", b.Ucell()[0]/set.NX(), b.Ucell()[1]/set.NX(), b.Ucell()[2]/set.NX(), b.Ucell()[3]/set.NY(), b.Ucell()[4]/set.NY(), b.Ucell()[5]/set.NY(), b.Ucell()[6]/set.NZ(), b.Ucell()[7]/set.NZ(), b.Ucell()[8]/set.NZ()); } file.Printf("#%s %s %s %s\n", Xdim.Label().c_str(), Ydim.Label().c_str(), Zdim.Label().c_str(), set.legend()); } std::string xyz_fmt; if (XcolPrecSet()) { TextFormat nfmt( XcolFmt(), XcolWidth(), XcolPrec() ); xyz_fmt = nfmt.Fmt() + " " + nfmt.Fmt() + " " + nfmt.Fmt() + " "; } else { TextFormat xfmt( XcolFmt(), set.NX(), Xdim.Min(), Xdim.Step(), 8, 3 ); TextFormat yfmt( XcolFmt(), set.NY(), Ydim.Min(), Ydim.Step(), 8, 3 ); TextFormat zfmt( XcolFmt(), set.NZ(), Zdim.Min(), Zdim.Step(), 8, 3 ); xyz_fmt = xfmt.Fmt() + " " + yfmt.Fmt() + " " + zfmt.Fmt() + " "; } if (sparse_) { for (pos[2] = 0; pos[2] < set.NZ(); ++pos[2]) { for (pos[1] = 0; pos[1] < set.NY(); ++pos[1]) { for (pos[0] = 0; pos[0] < set.NX(); ++pos[0]) { double val = set.GetElement(pos[0], pos[1], pos[2]); if (val > cut_) { Vec3 xyz = set.Bin().Corner(pos[0], pos[1], pos[2]); file.Printf( xyz_fmt.c_str(), xyz[0], xyz[1], xyz[2] ); set.WriteBuffer( file, pos ); file.Printf("\n"); } } } } } else { for (pos[2] = 0; pos[2] < set.NZ(); ++pos[2]) { for (pos[1] = 0; pos[1] < set.NY(); ++pos[1]) { for (pos[0] = 0; pos[0] < set.NX(); ++pos[0]) { Vec3 xyz = set.Bin().Corner(pos[0], pos[1], pos[2]); file.Printf( xyz_fmt.c_str(), xyz[0], xyz[1], xyz[2] ); set.WriteBuffer( file, pos ); file.Printf("\n"); } } } } return 0; }
// DataIO_Std::WriteSet2D() int DataIO_Std::WriteSet2D( DataSet const& setIn, CpptrajFile& file ) { if (setIn.Ndim() != 2) { mprinterr("Internal Error: DataSet %s in DataFile %s has %zu dimensions, expected 2.\n", setIn.legend(), file.Filename().full(), setIn.Ndim()); return 1; } DataSet_2D const& set = static_cast<DataSet_2D const&>( setIn ); int xcol_width = 8; int xcol_precision = 3; Dimension const& Xdim = static_cast<Dimension const&>(set.Dim(0)); Dimension const& Ydim = static_cast<Dimension const&>(set.Dim(1)); if (Xdim.Step() == 1.0) xcol_precision = 0; DataSet::SizeArray positions(2); TextFormat ycoord_fmt(XcolFmt()), xcoord_fmt(XcolFmt()); if (square2d_) { // Print XY values in a grid: // x0y0 x1y0 x2y0 // x0y1 x1y1 x2y1 // x0y2 x1y2 x2y2 // If file has header, top-left value will be '#<Xlabel>-<Ylabel>', // followed by X coordinate values. if (writeHeader_) { ycoord_fmt.SetCoordFormat( set.Nrows(), Ydim.Min(), Ydim.Step(), xcol_width, xcol_precision ); std::string header; if (Xdim.Label().empty() && Ydim.Label().empty()) header = "#Frame"; else header = "#" + Xdim.Label() + "-" + Ydim.Label(); WriteNameToBuffer( file, header, xcol_width, true ); xcoord_fmt.SetCoordFormat( set.Ncols(), Xdim.Min(), Xdim.Step(), set.Format().ColumnWidth(), xcol_precision ); for (size_t ix = 0; ix < set.Ncols(); ix++) file.Printf( xcoord_fmt.fmt(), set.Coord(0, ix) ); file.Printf("\n"); } for (positions[1] = 0; positions[1] < set.Nrows(); positions[1]++) { if (writeHeader_) file.Printf( ycoord_fmt.fmt(), set.Coord(1, positions[1]) ); for (positions[0] = 0; positions[0] < set.Ncols(); positions[0]++) set.WriteBuffer( file, positions ); file.Printf("\n"); } } else { // Print X Y Values // x y val(x,y) if (writeHeader_) file.Printf("#%s %s %s\n", Xdim.Label().c_str(), Ydim.Label().c_str(), set.legend()); if (XcolPrecSet()) { xcoord_fmt = TextFormat(XcolFmt(), XcolWidth(), XcolPrec()); ycoord_fmt = xcoord_fmt; } else { xcoord_fmt.SetCoordFormat( set.Ncols(), Xdim.Min(), Xdim.Step(), 8, 3 ); ycoord_fmt.SetCoordFormat( set.Nrows(), Ydim.Min(), Ydim.Step(), 8, 3 ); } std::string xy_fmt = xcoord_fmt.Fmt() + " " + ycoord_fmt.Fmt() + " "; for (positions[1] = 0; positions[1] < set.Nrows(); ++positions[1]) { for (positions[0] = 0; positions[0] < set.Ncols(); ++positions[0]) { file.Printf( xy_fmt.c_str(), set.Coord(0, positions[0]), set.Coord(1, positions[1]) ); set.WriteBuffer( file, positions ); file.Printf("\n"); } } } return 0; }
// Cluster_Kmeans::ClusterResults() void Cluster_Kmeans::ClusterResults(CpptrajFile& outfile) const { outfile.Printf("#Algorithm: Kmeans nclusters %i maxit %i\n", nclusters_, maxIt_); }
/** Write data at frame to CharBuffer. If no data for frame write 0.0. */ void DataSet_double::WriteBuffer(CpptrajFile &cbuffer, SizeArray const& frame) const { if (frame[0] >= Data_.size()) cbuffer.Printf(format_.fmt(), 0.0); else cbuffer.Printf(format_.fmt(), Data_[frame[0]]); }