Main Page | Modules | Data Structures | File List | Data Fields | Globals

mqparser/mqscriptfunc.h File Reference

A Documented file. More...

Go to the source code of this file.

MQScript Auxiliary Functions

double length (VariantArray array)
 Returns the length of an array.

double verbose (void)
 Toggles the verbose mode.

double Number_Save (double number, char *filename)
 Saves a number to a file.


MQScript Map2D Input/Output Functions

Map2D Map_Load (char *mqf_map)
 Loads and returns a data structure of a 2D-map (Map2D).

Map2D lmap (char *mqf_map)
 Short-hand version of Map_Load.

double Map_Save (Map2D map, char *mqf_map)
 Saves a 2D-map (Map2D) at the location filename.

double smap (Map2D map, char *mqf_map)
 Short-hand version of Map_Save.

double Map_Remove (char *mqf_map)
 Deletes a 2D-map file at the location mqf_map.

double rmmap (char *mqf_map)
 Short-hand version of Map_Remove.


MQScript Map2D Get Functions

double Map_GetArea (Map2D map)
 Calculates and returns the area of the 2D-map.

double Map_GetMedian (Map2D map, double mode)
 Calculates and returns the median of the abundance values for all the points in a 2D-map.

double Map_GetAvgAbsDevFromValue (Map2D map, double value, double mode)
 Calculates and returns the average absolute deviation from a given value for all the abundance values of the points in a 2D-map.

double Map_GetMedianAbsDevFromValue (Map2D map, double value, double mode)
 Calculates and returns the median absolute deviation from a given value for all the abundance values of the points in a 2D-map.

double Map_GetMean (Map2D map, double mode)
 Calculates and returns the mean of the abundance values for all the points in a 2D-map.

double Map_GetSD (Map2D map, double mean, double mode)
 Calculates and returns the standard deviation of the abundance values for all the points in a 2D-map.

VariantArray Map_GetSegmentArray (Map2D map)
 Returns all the segments of a segmented map in an array.

VariantArray Map_GetSegmentArrayContainingPeaks (Map2D map, VariantArray peaks, double dMZTolerance, double nScanTolerance)
 Returns only the segments of a segmented map that contain the peaks of interest.


MQScript Map2D Processing Functions

Map2D Map_ApplyFilter (Map2D map, char *filter, double dim)
 Applies a convolution filter to a 2D-map (Map2D).

Map2D filtermap (Map2D map, char *filter, double dim)
 Short-hand version of Map_ApplyFilter.

Map2D Map_Segment (Map2D map, double bin_size)
 Applies the watershed algorithm to a 2D-map (Map2D).

Map2D segmap (Map2D map, double bin_size)
 Short-hand version of Map_Segment.

Map2D Map_ApplyOpening (Map2D map, char *structel)
 Applies a morphological opening to the 2D-map (Map2D).

Map2D openmap (Map2D map, char *structel)
 Short-hand version of Map_ApplyOpening.

Map2D Map_ApplyClosing (Map2D map, char *structel)
 Applies a morphological closing to the 2D-map (Map2D).

Map2D closemap (Map2D map, char *structel)
 Short-hand version of Map_ApplyClosing.

Map2D Map_ApplyLinearInterpolation (Map2D map)
 Applies a linear interpolation for each mass chromatogram.

Map2D linintermap (Map2D map)
 Short-hand version of Map_ApplyLinearInterpolation.

Map2D Map_RemoveStreaks4 (Map2D map)
 Removes the streaks from a map.

Map2D rstreaks4 (Map2D map)
 Short-hand version of Map_RemoveStreaks4.

VariantArray Map_FindPeaks (Map2D source, char *structel, double abuthr)
 Finds peaks in a map.

VariantArray findpk (Map2D source, char *structel, double abuthr)
 Short-hand version of Map_FindPeaks.

Map2D Map_Extract (Map2D source, Map2D stencil, Segment segment)
 Extracts a segment map.

Map2D extractmap (Map2D source, Map2D stencil, Segment segment)
 Short-hand version of Map_Extract.


MQScript Peak Related Functions

VariantArray PeakArray_Load (char *mqf_pk)
 Loads a peak series.

VariantArray lpk (char *mqf_pk)
 Short-hand version of PeakArray_Load.

double PeakArray_Print (VariantArray array)
 Prints a peak series.

double PeakArray_Save (VariantArray array, char *mqf_pk)
 Saves a peak series.

double PeakArray_SaveWithMapValue (VariantArray array, Map2D map, char *mqf_pk)
 Saves a peak series with value in map.

double PeakArray_SavePhysQuantWithMapValue (VariantArray array, Map2D map, char *mqf_pk)
 Saves a peak series using their physical quantity values (time and m/z) and the intensity value.

double PeakArray_SaveWithBaseline (VariantArray array, Map2D map, char *mqf_pk)
 Saves a peak series with baseline values.


MQScript FPeak Related Functions

VariantArray FPeakArray_Load (char *mqf_fpk)
 Loads a fpeak series.

VariantArray lfpk (char *name)
 Short-hand version of FPeakArray_Load.

double FPeakArray_Print (VariantArray array)
 Prints an array of type FPeak2i to standard output.

VariantArray FPeakArray_GetPeaksNotInPeakGroup (VariantArray pFPeaksAll, PeakGroup peakgroup)
 Returns the subset of fpeaks from a set of fpeaks, that do not belong to a particular peakgroup.

VariantArray FPeakArray_GetPeakArray (VariantArray array)
 Converts an array of type FPeak2i to an array of type Peak2i.

Map2D FPeakArray_CreateMapFrom (VariantArray array, Map2D template_map)
 Creates a map from an array of peaks of type FPeak2i.


MQScript Segment Related Functions

Segment Segment_FindAndFitPeaks (Segment segment, Map2D map, char *structel, double abuthr, char *curve)
 Finding and fitting peaks in a 2D-map.

Map2D Segment_GetMap (Segment segment, Map2D source, Map2D mask)
 Extracts a segment map.

double Segment_Save (Segment segment, char *mqf_fpk)
 Saves a Segment in a file mqf_fpk.

double Segment_Print (Segment segment)
 Prints a Segment to standard output.

VariantArray Segment_GetFPeaks (Segment segment)
 Gets the FPeak2i array from a Segment.

VariantArray SegmentArray_Load (char *mqf_fpk)
 Loads an array of type Segment.

double FPeakFile_WriteHeader (char *szFilename, Map2D source, char *curve, double nsegments)
 Saves the header of a fpeaks file.

double FPeakFile_CalcStatAndSave (char *mqf_fpk)
 Calculates and saves the statistics of attributes of the fpeaks.


MQScript PeakGroup Related Functions

VariantArray PeakGroupArray_Load (char *mqf_pgr)
 Loads a peakgroup array from a pgroups file.

VariantArray lpg (char *mqf_pgr)
 Short-hand version of PeakGroupArray_Load.

double PeakGroup_Print (PeakGroup peakgroup)
 Prints a PeakGroup in standard out.

double PeakGroup_Save (PeakGroup peakgroup, char *filename, double type)
 Saves a PeakGroup in a file.

VariantArray PeakGroup_GetFPeaks (PeakGroup peakgroup)
 Gets the FPeak2i array from a PeakGroup.

VariantArray getfpk (PeakGroup peakgroup)
 Short-hand version of PeakGroup_GetFPeaks.

PeakGroup PeakGroup_Create (VariantArray array)
 Creates a peakgroup.

PeakGroup mkpg (VariantArray array)
 Short-hand version of PeakGroup_Create.

VariantArray PeakGroupArray_GetFPeaks (VariantArray pg_array)
 Gets the FPeak2i array from a PeakGroup array.

VariantArray vgetfpk (VariantArray pg_array)
 Short-hand version of PeakGroupArray_GetFPeaks.

Map2D getpgmap (PeakGroup peakgroup, Map2D source, Map2D stencil)
 Returns the peakgroup map of a peakgroup.

VariantArray getpgmappeaks (PeakGroup peakgroup, char *mqf_fpk)
 Returns the fpeaks that are found in a peakgroup map.

double PeakGroupFile_WriteHeader (char *mqf_fpk, char *type, char *curve, double npg, double npeaks)
 Saves the header of a pgroups file.


MQScript FPeak Clustering Related Functions

double FPeakFile_ClusterAndSave (char *mqf_fpk_in, char *mgf_pgr_out, double scan_thr, double mz_thr, double abund_thr, double nLinkageType)
 Clusters fpeaks (FPeak2i) into co-eluting peakgroups (PeakGroup).

double PeakGroupFile_ClusterAndSave (char *mgf_pgr_in, char *mgf_pgr_out, double scan_thr, double mz_thr, double abund_thr, double nLinkageType)
 Clusters fpeaks (FPeak2i) iteratively in each PeakGroup of a PeakGroup array.


MQScript SegmentHandler Related Functions

SegmentHandler SegmentHandler_Init (char *peakgroupfilename)
SegmentHandler SegmentHandler_Set (SegmentHandler seghand, PeakGroup peakgroup)
Map2D SegmentHandler_GetMap (SegmentHandler seghand, PeakGroup peakgroup)
VariantArray SegmentHandler_GetMapPeaks (SegmentHandler seghand, PeakGroup peakgroup)

MQScript PeakGroup Refining/Deconvolving Functions

PeakGroup PeakGroup_Refine7 (PeakGroup peakgroup, Map2D PGMap, VariantArray pPGOutFPeaks, char *struct_elem, double abu_thr, char *curve, double scan_window, double mz_overlap_thr)
 Deconvolves overlapping fpeaks in a PeakGroup.

PeakGroup refine7 (PeakGroup peakgroup, Map2D PGMap, VariantArray pPGOutFPeaks, char *struct_elem, double abu_thr, char *curve, double scan_window, double mz_overlap_thr)
 Short-hand version of PeakGroup_Refine7.

PeakGroup PeakGroup_Deconv3 (PeakGroup peakgroup, Map2D PGMap, VariantArray pPGOutFPeaks, char *curve)
 Deconvolves the PeakGroup into IsotopicClusters using fitting.

PeakGroup deconv3 (PeakGroup peakgroup, Map2D PGMap, VariantArray pPGOutFPeaks, char *curve)
 Short-hand version of PeakGroup_Deconv3.

PeakGroup PeakGroup_Deconv4 (PeakGroup peakgroup)
 Deconvolves the PeakGroup into IsotopicClusters without using fitting.

PeakGroup deconv4 (PeakGroup peakgroup)
 Short-hand version of PeakGroup_Deconv4.


MQScript Manual/Visual Fitting Related Functions

double fixmz (void)
double fixz (void)
double fixc (void)
double fixsmz (void)
double setsmz (double index, double value)
double clrfitter (void)
double commonsmz (void)
VariantArray fit (double init_sigma_rt, double init_sigma_mz)
VariantArray fitsmz (double init_sigma_rt, double smz_begin, double smz_end, double delta_sigma_mz)
VariantArray fitpep (double charge, double carbons)
double sfpk (VariantArray array, char *mqf)

MQScript Visualization Related Functions

double plotmap (Map2D map)
 Plots a variable of type Map2D.

double overlaymap (Map2D map)
 Overlays the segment boundaries of a labelled map on a grayscale map.

double plotpk (VariantArray array)
 Plots an array of peaks, either of type Peak2i or FPeak2i.

double plotpg (VariantArray array)
 Plots an array of peakgroups.

double plotseqsum (char *filename, char *which)
 Plots the sequest hits summarized in a seqsum file.

double movepk (double handle, double dx)
 Moves the rendering of the peaks stored at the given handle position.

double plotmqsqmap (char *filename)
 Plots the contents of an mqsqmap file.

double plotms2 (char *ms2_exp_name)
 Plots all the ms/ms (sequencing events) on the map.

double coversegments (VariantArray peak_array, double mz_threshold, double scan_theshold)
 Covers the segments that overlap with any of the rectangles surrounding the array of peaks.

double coversegmentsonedge (void)
 Covers the segments that touch the top and bottom edges of a Map2D.

double showlabels (void)
 Plots the label number of each segment.

double axes (void)
 Toggles the values of the axes from sampling coordinates to physical quantities.

double grid (void)
 Toggles the view of the grid.

double zoomin (double x1, double x2, double y1, double y2)
 Zooms in the required rectangle.

double setgraylevel (char *level)
 Sets the graylevel of the map.

double goscan (double scan)
 Displays the mass spectrum of that scan.

double show (double handle, char *field)
 Shows.

double hide (double handle)
 Hides.

double clear (double handle)
 Clears.


MQScript Other Functions

Map2D extractseam (Map2D smap1, Map2D smap2, Map2D lmap1, Map2D lmap2)
double tilefpeakfile (char *mqf_fpk, char *mqf_pgf)
char * print (char *string)
char * int2str (double number, double width)
double setfilter (char *filter)
double createMS2table (char *expmntname)
double printZoomBoundaries (void)
double auxExtractMapShown (char *filename)
double auxExtractMap (Map2D map, double x1, double x2, double y1, double y2, char *filename)
double PrintFPeaksArrayInfo (VariantArray fpeaks, double parameter)
double PeakArray_SaveDTAPositionsWithSequence (char *szSequestSummaryName, char *szSExperimentName, char *szQExperimentName, char *szOutputfile)
VariantArray PeakArray_GetDTAPositions (char *szSequestSummaryName, char *szSExperimentName, char *szQExperimentName)
double MS2GuidedQuantitation (Map2D map, double bin_width)


Detailed Description

A Documented file.

Author:
Kyriacos C. Leptos leptos@fas.harvard.edu
Date:
2004-10-05
Details.

Function Documentation

double axes void   ) 
 

Toggles the values of the axes from sampling coordinates to physical quantities.

This function when called requires that window that has a map already plotted in it

double clear double  handle  ) 
 

Clears.

Parameters:
handle 

double coversegments VariantArray  peak_array,
double  mz_height,
double  scan_theshold
 

Covers the segments that overlap with any of the rectangles surrounding the array of peaks.

This function when called requires that window that has a grayscale map already plotted in it

Parameters:
peak_array the array of type Peak2i
mz_threshold half of the height of the window surrounding each peak along the spectrometer dimension
scan_threshold half of the width of the window surrounding each peak along the chromatography dimension

double coversegmentsonedge void   ) 
 

Covers the segments that touch the top and bottom edges of a Map2D.

This function when called requires that window that has a grayscale map already plotted in it

Map2D FPeakArray_CreateMapFrom VariantArray  array,
Map2D  template_map
 

Creates a map from an array of peaks of type FPeak2i.

Parameters:
array the array of type FPeak2i
template_map this is a data structure of type Map2D that holds the header infromation for the newly constructed map

VariantArray FPeakArray_GetPeakArray VariantArray  array  ) 
 

Converts an array of type FPeak2i to an array of type Peak2i.

Returns an array of type Peak2i that was generated by stripping all the fitting parameters from each peak of type FPeak2i from the variable array

Parameters:
array the array of type FPeak2i to be converted

VariantArray FPeakArray_GetPeaksNotInPeakGroup VariantArray  pFPeaksAll,
PeakGroup  peakgroup
 

Returns the subset of fpeaks from a set of fpeaks, that do not belong to a particular peakgroup.

This function returns an array of type FPeak2i as a subset from the fpeak array defined by the variable pFPeaksAll, that do not belong to the PeakGroup varaible peakgroup

Parameters:
pFPeaksAll the array of type FPeak2i
peakgroup the PeakGroup under investigation

VariantArray FPeakArray_Load char *  mqf_fpk  ) 
 

Loads a fpeak series.

Loads and returns an array of fpeaks (FPeak2i)

Parameters:
mqf_fpk the name of the fpeaks file to be loaded The name uses the MapQuant naming convention

double FPeakArray_Print VariantArray  array  ) 
 

Prints an array of type FPeak2i to standard output.

Parameters:
array an array of fpeaks (FPeak2i) to be printed

double FPeakFile_CalcStatAndSave char *  mqf_fpk  ) 
 

Calculates and saves the statistics of attributes of the fpeaks.

This function loads an array of FPeak2i from a fpeaks file given as mqf_fpk. It then calculates and saves the statistics of attributes of the fpeaks into a file with the same name as mqf_fpk but replacing the extention to metrics

Parameters:
mqf_fpk the name of the fpeaks file to be loaded The filename uses the MapQuant naming convention

double FPeakFile_ClusterAndSave char *  mqf_fpk_in,
char *  mgf_pgr_out,
double  scan_thr,
double  mz_thr,
double  abund_thr,
double  nLinkageType
 

Clusters fpeaks (FPeak2i) into co-eluting peakgroups (PeakGroup).

This function reads the fpeaks from the file szInFilename, clusters them into peak groups and then outputs them into the file szOutFilename. The arguments scan_thr is used to set the scan tolerance of the co-elution restriction mentioned above. Likewise the argument mz_thr is used to set the m/z tolerance for the maximum allowed distance between two isotopic peaks.

Parameters:
mqf_fpk_in the fpeaks file from where the fpeaks to be clustered are loaded The filename uses the MapQuant naming convention
mqf_pgr_out the pgroups file to where the peakgroups are to be saved The filename uses the MapQuant naming convention
scan_thr the maximum distance along the chromatography dimension (in scan units) that two neighboring peaks in a cluster can have
mz_thr the maximum distance along the m/z (in m/z units) that two neighboring peaks in a cluster can have
abund_thr the minimum abundance that a peak needs to have in order to be condifered in the clustering
nLinkageType the type of clustering used CLUSTER_SINGLE_LINK or CLUSTER_HYBRID_LINK

double FPeakFile_WriteHeader char *  mqf_fpk,
Map2D  source,
char *  curve,
double  nsegments
 

Saves the header of a fpeaks file.

Parameters:
mqf_fpk the filename of the fpeaks file The filename uses the MapQuant naming convention
source the grayscale map used for peak finding and fitting
curve the curve used to fit the peaks
nsegments number of segments to be saved in an fpeaks file

Map2D getpgmap PeakGroup  peakgroup,
Map2D  source,
Map2D  stencil
 

Returns the peakgroup map of a peakgroup.

The peakgroup map of peakgroup is defined as the union of all the segment maps that contain the fpeaks found in the peakgroup itself. Segment maps are extracted as in the function Map_Extract using a Map2D as a source and a Map2D as stencil.

Parameters:
peakgroup the name of the file to be loaded
source the grayscale map used for extraction
stenicl the stencil map used that defines the points that belong to the segment

VariantArray getpgmappeaks PeakGroup  peakgroup,
char *  mqf_fpk
 

Returns the fpeaks that are found in a peakgroup map.

This function loads all the fpeaks (FPeak2i) from the file mqf_fpk and returns the all fpeaks that are found in the segments that comprise the corresponding peakgroup 2D-map defined by peakgroup. The peakgroup map of PeakGroup is defined as the union of all the segment maps that contain the fpeaks found in the peakgroup itself.

Parameters:
peakgroup the name of the PeakGroup udner investigation
mqf_fpk the name of an fpeaks file using the MapQuant naming convention

double goscan double  scan  ) 
 

Displays the mass spectrum of that scan.

This function when called requires that window that has a grayscale map already plotted in it

Parameters:
scan the scan to be displayed

double grid void   ) 
 

Toggles the view of the grid.

This function when called requires that window that has a map already plotted in it

double hide double  handle  ) 
 

Hides.

Parameters:
handle 
field 

double length VariantArray  array  ) 
 

Returns the length of an array.

Returns the length of a VariantArray

Parameters:
array the array whose length is going to be returned

Map2D Map_ApplyClosing Map2D  map,
char *  structel
 

Applies a morphological closing to the 2D-map (Map2D).

Morphological closing is an image operation that is used to fill in the holes for groups of points (pixels) whose connectivity is defined by the structuring element given in structel

Parameters:
map the 2D-map being processed
structel a string definining a particular structuring element

Map2D Map_ApplyFilter Map2D  map,
char *  filter,
double  dim
 

Applies a convolution filter to a 2D-map (Map2D).

Saves a 2D-map (Map2D) as the name of variable filename indicates

Parameters:
map the 2D-map being processed
filter A string representing the convolution kernel for the noise filtering The value of filter can either be "BC_n" or "GS_n_m" where n and m are integers. BC refers to a box-car filter (otherwise called moving-average filter) of n points and GS refers to a Gaussian filter of standard deviation n.m. Savitzky-Golay filters can also be applied in the form of "SG_m_n_M" where m is points to the left, n is points to the right and M is the moment.
dim The dimension to which the filter is going to be applied. It takes the values RT_DIMENSION for the chromatography dimension and MZ_DIMENSION referring to the dimension of the 2-D map that the filter will be applied.

Map2D Map_ApplyLinearInterpolation Map2D  map  ) 
 

Applies a linear interpolation for each mass chromatogram.

Applies a linear interpolation for each mass chromatogram for scans that are not MS, i.e. dependent scans

Parameters:
map the 2D-map being processed

Map2D Map_ApplyOpening Map2D  map,
char *  structel
 

Applies a morphological opening to the 2D-map (Map2D).

Morphological opening is an image operation that is used to select for groups of points (pixels) whose connectivity is defined by the structuring element given in structel

Parameters:
map the 2D-map being processed
structel a string definining a particular structuring element

Map2D Map_Extract Map2D  source,
Map2D  stencil,
Segment  segment
 

Extracts a segment map.

A function that excises and returns a sub Map2D from a source Map2D whose borders are defined by the rectangle that circumscribes the Segment in query. The function uses the Map2D mask as a stencil.

Parameters:
source the grayscale map used for extraction
stencil the stencil map used that defines the points that belong to the segment
segment a segment data sctructure containing the borders of the 2D-map that circumscribes the actual segment and used as the primary information to excise the segment 2D-map

VariantArray Map_FindPeaks Map2D  source,
char *  structel,
double  abuthr
 

Finds peaks in a map.

A function that finds and returns an array of Peak2i. It uses the structuring element structel and the noise threshold abuthr

Parameters:
source the grayscale map used for peak finding
structel the string defining a structuring element for the peak-finding algorithm
abuthr the abundance threshold used for peak finding

double Map_GetArea Map2D  map  ) 
 

Calculates and returns the area of the 2D-map.

Calculates the area of the map using the SCAN_SIZE and BIN_SIZE values

Parameters:
map a Map2D data-structure

double Map_GetAvgAbsDevFromValue Map2D  map,
double  value,
double  mode
 

Calculates and returns the average absolute deviation from a given value for all the abundance values of the points in a 2D-map.

Calculates the average absolute deviation from a given value for the values int the abundance matrix of the 2D-map, either by including zeros values or excluding them.

Parameters:
map a Map2D data-structure

double Map_GetMean Map2D  map,
double  mode
 

Calculates and returns the mean of the abundance values for all the points in a 2D-map.

Calculates the statistical mean of the abundance matrix of the 2D-map, either by including zeros values or excluding them.

Parameters:
map a Map2D data-structure

double Map_GetMedian Map2D  map,
double  mode
 

Calculates and returns the median of the abundance values for all the points in a 2D-map.

Calculates the statistical median of the abundance matrix of the 2D-map, either by including zeros values or excluding them.

Parameters:
map a Map2D data-structure

double Map_GetMedianAbsDevFromValue Map2D  map,
double  value,
double  mode
 

Calculates and returns the median absolute deviation from a given value for all the abundance values of the points in a 2D-map.

Calculates the median absolute deviation from a given value for the values int the abundance matrix of the 2D-map, either by including zeros values or excluding them.

Parameters:
map a Map2D data-structure

double Map_GetSD Map2D  map,
double  mean,
double  mode
 

Calculates and returns the standard deviation of the abundance values for all the points in a 2D-map.

Calculates the standard deviation of the abundance matrix of the 2D-map, either by including zeros values or excluding them.

Parameters:
map a Map2D data-structure

VariantArray Map_GetSegmentArray Map2D  map  ) 
 

Returns all the segments of a segmented map in an array.

When a 2D-map is segmented using watershed segmentation, information of the boundaries of each segment are stored internally in an array of type Segment

Parameters:
map a Map2D data-structure

VariantArray Map_GetSegmentArrayContainingPeaks Map2D  map,
VariantArray  peaks,
double  dMZTolerance,
double  nScanTolerance
 

Returns only the segments of a segmented map that contain the peaks of interest.

Similar to the function Map_GetSegmentArray but restricting it to the segments that are found within the areas defined by the rectangles defined by the parameters dMZTolerance and nScanTolerance

Parameters:
map a Map2D data-structure
peaks an array of Peak2i that define the centers of the rectangles of interest
dMZTolerance half the size of the length of the rectangle of interest surrounding each peak
nScanTolerance half the size of the width of the rectangle of interest surrounding each peak

Map2D Map_Load char *  mqf_map  ) 
 

Loads and returns a data structure of a 2D-map (Map2D).

Loads and returns a data structure of a 2D-map (Map2D). The filename should can also be in the MQScript convention format

Parameters:
mqf_map the filename of the 2D-map that is being loaded using the MapQuant naming convention

double Map_Remove char *  mqf_map  ) 
 

Deletes a 2D-map file at the location mqf_map.

Saves a 2D-map (Map2D) as the name of variable filename indicates

Parameters:
mqf_map the location where the 2D-map file is stored using the MapQuant naming convention

Map2D Map_RemoveStreaks4 Map2D  map  ) 
 

Removes the streaks from a map.

Removes streaks from a Map2D structure by subtracting the median value of each mass chromatogram from all the points in the mass chromatogram

Parameters:
map the 2D-map being processed

double Map_Save Map2D  map,
char *  mqf_map
 

Saves a 2D-map (Map2D) at the location filename.

Saves a 2D-map (Map2D) as the name of variable filename indicates

Parameters:
map the 2D-map being saved
mqf_map the filename of the 2D-map that is being saved

Map2D Map_Segment Map2D  map,
double  bin_size
 

Applies the watershed algorithm to a 2D-map (Map2D).

This function takes the 2-D map map and segments it using the watershed algorithm. The watershed algorithm is a segmentation algorithm that divides the map into regions where the signal is concentrated thus rendering the analysis of the map easier. The parameter bin_size is used to denote the ion abundance stepping size used by the flooding procedure of the watershed algorithm. The function returns a 2-D map that can be used a mask, since each of its data points holds the segment number in which it belongs to. It returns an array of type Segment.

Parameters:
map the 2D-map being processed
bin_size the step (in intensity units) used to flood the map in order to build the watersheds (boundaries)

double movepk double  handle,
double  dx
 

Moves the rendering of the peaks stored at the given handle position.

This function when called requires that window that has a grayscale map already plotted in it

Parameters:
handle 
dx 

double Number_Save double  number,
char *  filename
 

Saves a number to a file.

Saves a number (in append mode) to the file named filename.

Parameters:
number the number to be saved
filename the filename where the number is being saved

double overlaymap Map2D  map  ) 
 

Overlays the segment boundaries of a labelled map on a grayscale map.

This function when called requires that window that has a grayscale map already plotted in it

Parameters:
map a Map2D data-structure that is loaded either from the command line or by using the lmap function

VariantArray PeakArray_Load char *  mqf_fpk  ) 
 

Loads a peak series.

Loads a plain peak series (Peak2i) into an array data structure

Parameters:
mqf_fpk the name of the file to be loaded

double PeakArray_Print VariantArray  array  ) 
 

Prints a peak series.

Prints an array of type Peak2i in standard out

Parameters:
array the array if type Peak2i to be printed

double PeakArray_Save VariantArray  array,
char *  mqf_pk
 

Saves a peak series.

Saves a peak series in a file

Parameters:
array an array of type Peak2i to be saved
mqf_pk the filename where the data are saved. The filename uses the MapQuant naming convention

double PeakArray_SavePhysQuantWithMapValue VariantArray  array,
Map2D  map,
char *  mq_pk
 

Saves a peak series using their physical quantity values (time and m/z) and the intensity value.

Saves a peak series as 3-column table in a file. The columns in the file are the following: 1. retention time of the peak 2. m/z of the peak 3. abundance in map of that (rt, m/z) point

Parameters:
array an array of type Peak2i to be saved
map the Map2D structure from where the intensity values are taken
mqf_pk the filename where the data are saved The filename uses the MapQuant naming convention

double PeakArray_SaveWithBaseline VariantArray  array,
Map2D  map,
char *  mqf_pk
 

Saves a peak series with baseline values.

Saves a peak series as 4-column table in a file. The columns in the file are the following: 1. retention time of the peak 2. m/z of the peak 3. scan position defining the left shoulder of the peak 4. scan position defining the right shoulder of the peak

Parameters:
array an array of type Peak2i to be saved
map the Map2D structure from where the baseline position values are taken
mqf_pk the filename where the data are saved The filename uses the MapQuant naming convention

double PeakArray_SaveWithMapValue VariantArray  array,
Map2D  map,
char *  mqf_pk
 

Saves a peak series with value in map.

Saves a peak series as 3-column table in a file. The columns in the file are the following: 1. scan position of the peak 2. m/z bin position of the peak 3. abundance in map at that (scan, m/z bin) point

Parameters:
array an array of type Peak2i to be saved
map the Map2D structure from where the intensity values are taken
mqf_pk the filename where the data are saved The filename uses the MapQuant naming convention

PeakGroup PeakGroup_Create VariantArray  array  ) 
 

Creates a peakgroup.

Parameters:
array 

PeakGroup PeakGroup_Deconv3 PeakGroup  peakgroup,
Map2D  PGMap,
VariantArray  pPGOutFPeaks,
char *  curve
 

Deconvolves the PeakGroup into IsotopicClusters using fitting.

Parameters:
peakgroup the PeakGroup to be deconvolved
PGMap the peakgroup map used for the deconvolution
pPGOutFPeaks the fpeaks that do not belong in the peakgroup
curve the curve to be used to fit the peaks

PeakGroup PeakGroup_Deconv4 PeakGroup  peakgroup  ) 
 

Deconvolves the PeakGroup into IsotopicClusters without using fitting.

Parameters:
peakgroup the PeakGroup to be deconvolved

VariantArray PeakGroup_GetFPeaks PeakGroup  peakgroup  ) 
 

Gets the FPeak2i array from a PeakGroup.

Parameters:
peakgroup the PeakGroup from which the FPeak2i array is taken

double PeakGroup_Print PeakGroup  peakgroup  ) 
 

Prints a PeakGroup in standard out.

Parameters:
peakgroup the PeakGroup to be printed

PeakGroup PeakGroup_Refine7 PeakGroup  peakgroup,
Map2D  PGMap,
VariantArray  pPGOutFPeaks,
char *  struct_elem,
double  abu_thr,
char *  curve,
double  scan_window,
double  mz_overlap_thr
 

Deconvolves overlapping fpeaks in a PeakGroup.

Parameters:
peakgroup the peakgroup to be refined
PGMap the peakgroup map used in the deconvolution
pPGOutFPeaks the fpeaks that do not belong in the peakgroup
struct_elem the string defining a structuring element for the peak-finding algorithm
curve the curve to be used to fit the peaks
scan_window ????
mz_overlap_thr the m/z threshold that defines which

double PeakGroup_Save PeakGroup  peakgroup,
char *  mqf_pgr,
double  type
 

Saves a PeakGroup in a file.

Parameters:
peakgroup the PeakGroup to be saved
mqf_pgr the filename where the PeakGroup is to be saved The filename uses the MapQuant naming convention
type the values it can take is WITH_CURVE or WITHOUT_CURVE

VariantArray PeakGroupArray_GetFPeaks VariantArray  pg_array  ) 
 

Gets the FPeak2i array from a PeakGroup array.

Parameters:
pg_array an array of type PeakGroup

VariantArray PeakGroupArray_Load char *  mqf_pgr  ) 
 

Loads a peakgroup array from a pgroups file.

Parameters:
mqf_pgr the name of the pgroups file to be loaded

double PeakGroupFile_ClusterAndSave char *  mgf_pgr_in,
char *  mgf_pgr_out,
double  scan_thr,
double  mz_thr,
double  abund_thr,
double  nLinkageType
 

Clusters fpeaks (FPeak2i) iteratively in each PeakGroup of a PeakGroup array.

Parameters:
mqf_pgr_in the pgroups file from where the peakgroups whose fpeaks are going to be clustered are loaded The filename uses the MapQuant naming convention
mqf_pgr_out the pgroups file to where the peakgroups are to be saved The filename uses the MapQuant naming convention
scan_thr the maximum distance along the chromatography dimension (in scan units) that two neighboring peaks in a cluster can have
mz_thr the maximum distance along the m/z (in m/z units) that two neighboring peaks in a cluster can have
abund_thr the minimum abundance that a peak needs to have in order to be condifered in the clustering
nLinkageType the type of clustering used CLUSTER_SINGLE_LINK or CLUSTER_HYBRID_LINK

double PeakGroupFile_WriteHeader char *  mqf_fpk,
char *  type,
char *  curve,
double  npg,
double  npeaks
 

Saves the header of a pgroups file.

Parameters:
mqf_fpk the name of the file to be saved
type (NOT USED, ALWAYS USE "d3")
curve the curve to be used to fit the peaks
npg the number of peakgroups to be saved in the file
npeaks (NOT_USED, ALWAYS USE -1)

double plotmap Map2D  map  ) 
 

Plots a variable of type Map2D.

When this function if called from the mqshell a new window pops-up showing the map in grayscale

Parameters:
map a Map2D data-structure that is loaded either from the command line or by using the lmap function

double plotmqsqmap char *  filename  ) 
 

Plots the contents of an mqsqmap file.

This function when called requires that window that has a grayscale map already plotted in it

Parameters:
map a Map2D data-structure that is loaded either from the command line or by using the lmap function

double plotms2 char *  ms2_exp_name  ) 
 

Plots all the ms/ms (sequencing events) on the map.

This function when called requires that window that has a grayscale map already plotted in it

Parameters:
ms2_exp_name the name of an expreriment that contains

double plotpg VariantArray  array  ) 
 

Plots an array of peakgroups.

This function when called requires that window that has a grayscale map already plotted in it

Parameters:
array the

double plotpk VariantArray  array  ) 
 

Plots an array of peaks, either of type Peak2i or FPeak2i.

This function when called requires that window that has a grayscale map already plotted in it

Parameters:
array the

double plotseqsum char *  filename,
char *  which
 

Plots the sequest hits summarized in a seqsum file.

This function when called requires that window that has a grayscale map already plotted in it

Parameters:
filename the full path, name included, of the seqsum file
which this string acts as a filter, for the time being always use the string "all"

Segment Segment_FindAndFitPeaks Segment  segment,
Map2D  map,
char *  structel,
double  abuthr,
char *  curve
 

Finding and fitting peaks in a 2D-map.

Peak-finding
This function finds and fits peaks in the 2-D map map, and returns them as a variable of type Segment. The arguments used by the peak-finding part of the algorithm are structel, the structuring element and abuthr, the abundance cut off.

Peak-fitting
The argument used for the peak-fitting part of the algorithm is and it denotes the type of curve that the peak should be fitted to.

Parameters:
segment the segment to which the fpeaks found are going to be attached
map the gray-scale map used for peak finding and fitting
structel the string defining a structuring element for the peak-finding algorithm. It can take values such as "N9x3E" as shown here
abuthr the abundance threshold value which needs to be met by all the points in the structuring element and used for peak finding
curve a string describing the curve to be used to fit the peaks. It can take values as "NR_GAUSSIOID" (Gaussiod using Numerical Recipies algorithm)

\[f(m, r; A, r_o, m_o, \sigma_m, \sigma_r) = \frac{A}{2\pi\sigma_m\sigma_r}e^{\frac{(r-r_o)^2}{2\sigma_r^2}}e^{\frac{(m-m_o)^2}{2\sigma_m^2}}\]

or "NR_EM_GAUSSIOID" (Exponentailly Modified Gaussioid again using Numerical Recipies algorithm). Thes second curve models the trailing observed in chromatographic peaks, however this is only the beta version and has not been tested yet.

VariantArray Segment_GetFPeaks Segment  segment  ) 
 

Gets the FPeak2i array from a Segment.

Parameters:
segment a Segment from which the FPeak2i array are taken

Map2D Segment_GetMap Segment  segment,
Map2D  source,
Map2D  mask
 

Extracts a segment map.

Same as Map_Extract

Parameters:
segment a segment data sctructure containing the borders of the 2D-map that circumscribes the actual segment and used as the primary information to excise the segment 2D-map
source the grayscale map used for extraction
mask the stencil map used that defines the points that belong to the segment

double Segment_Print Segment  segment  ) 
 

Prints a Segment to standard output.

Parameters:
segment the Segment to be printed

double Segment_Save Segment  segment,
char *  mqf_fpk
 

Saves a Segment in a file mqf_fpk.

Parameters:
segment the Segment to be saved
mqf_fpk the filename of the fpeaks file where the Segment to be saved The filename uses the MapQuant naming convention

VariantArray SegmentArray_Load char *  mqf_fpk  ) 
 

Loads an array of type Segment.

Parameters:
mqf_fpk the filename of the fpeaks file from where an array of type Segment is to be loaded The filename uses the MapQuant naming convention

double setgraylevel char *  level  ) 
 

Sets the graylevel of the map.

This function when called requires that window that has a grayscale map already plotted in it

Parameters:
level a string that can take the values "log" or "linear"

double show double  handle,
char *  field
 

Shows.

Parameters:
handle 
field 

double showlabels void   ) 
 

Plots the label number of each segment.

This function when called requires that window that has a grayscale map with an stencil map overlayed on top of it already plotted in it

double verbose void   ) 
 

Toggles the verbose mode.

If called the verbose mode will be turned off or on according to the current state

double zoomin double  x1,
double  x2,
double  y1,
double  y2
 

Zooms in the required rectangle.

This function when called requires that window that has a grayscale map already plotted in it The filename uses the MapQuant naming convention

Parameters:
x1 left border (scan units)
x2 right border (scan units)
y1 lower border (m/z bin units)
y2 upper border (m/z bin units)


Generated on Sun Feb 13 01:06:02 2005 for MapQuant by doxygen 1.3.7