Structures |
---|
Name: Libnucnet__Species Description: Libnucnet__Species is a structure that stores nuclear data for a particular nuclear species. Routines act on the structure to retrieve data, update data, or use data to compute functions of the nuclear data, such as the nuclear partition function. The contents of the structure are not made public by the API. |
Name: Libnucnet__Nuc Description: Libnucnet__Nuc is a structure that stores the collection of nuclear species. Routines act on the structure to retrieve species or add or remove them. The contents of the structure are not made public by the API. |
Name: Libnucnet__NucView Description: Libnucnet__NucView is a structure that stores a view of Libnucnet__Nuc structure. A view is a subset of the parent Libnucnet__Nuc with the nuclei chosen by an XPath expression. The view structure contains a Libnucnet__Nuc pointer, which may be passed into all API routines that take such structures as input. It is important to note that the Libnucnet__Nuc member of a view does not own the nuclide, so modifying the data in the view modifies the data in the parent structure. |
User-Supplied Routines |
---|
Name: Libnucnet__Species__compare_function() Description: User-supplied routine to be applied during a sort of the species in the species collection. Syntax: int Libnucnet__Species__compare_function( const Libnucnet__Species *p_species1, const Libnucnet__Species *p_species2 );Input:
User's routine must return -1 if species 1 is less than species 2, 0 if the two species are equal, and 1 if species 1 is greater than species 2. If this routine is not supplied through Libnucnet__Nuc__setSpeciesCompareFunction, the default function, which sorts in increasing Z and A, is used. |
Name: Libnucnet__Species__iterateFunction() Description: User-supplied routine to be applied during an iteration over the species in the species collection. Syntax: int Libnucnet__Species__iterateFunction( Libnucnet__Species *self, void *p_user_data );Input:
User's routine must return 1 to continue or 0 to stop. |
Name: Libnucnet__Species__nseCorrectionFactorFunction() Description: Optional user-supplied routine to compute corrections to the NSE factor for a species. Syntax: double Libnucnet__Species__nseCorrectionFactorFunction( Libnucnet__Species *self, double d_t9, double d_rho, double d_ye, void *p_user_data );Input:
User's routine must return a double giving the NSE correction factor for the input temperature, density, and electron-to-baryon ratio. |
Routines |
---|
Name: Libnucnet__NucView__free() Description: Free the memory allocated for a view of a Libnucnet nuclear species collection. Syntax: void Libnucnet__NucView__free( Libnucnet__NucView *self );Input:
Upon successful return, the memory for the view has been freed. Example: Free up memory used for Libnucnet__NucView *p_my_view:Libnucnet__NucView__free( p_my_view ); |
Name: Libnucnet__NucView__getNuc() Description: Get the nuclear collection in a view. Syntax: Libnucnet__Nuc * Libnucnet__NucView__getNuc( Libnucnet__NucView *self );Input:
Routine returns the nuclear collection from a Libnucnet__NucView. If the input view is invalid, error handling is invoked. Example: Get the collection from Libnucnet__NucView *p_my_view:p_nuc = Libnucnet__NucView__getNuc( p_my_view ); |
Name: Libnucnet__NucView__new() Description: Create a new Libnucnet__NucView structure. Syntax: Libnucnet__NucView * Libnucnet__NucView__new( const Libnucnet__Nuc * p_nuc, const char * s_nuc_xpath );Input:
The routine returns a pointer to a new Libnucnet__NucView structure, that is, a view of a parent nuclide collection. If it is not possible to allocate memory for the new structure, Libnucnet error handling is invoked. Example: Create a view of the neon isotopes in Libnucnet__Nuc structure p_my_nuclei:p_view = Libnucnet__NucView__new( p_my_nuclei, "[z = 10]" ); |
Name: Libnucnet__Nuc__addSpecies() Description: Add a species to a species collection. Syntax: int Libnucnet__Nuc__addSpecies( Libnucnet__Nuc *self, Libnucnet__Species *p_species );Input:
Routine returns 1 for successful addition and 0 for failure. On successful return, the routine updates the Libnucnet__Nuc structure to include the new species. Example: Add p_species to the species collection p_my_nuclei:if( Libnucnet__Nuc__addSpecies( p_my_nuclei, p_species ) ) printf( "Addition successful.\n" ); |
Name: Libnucnet__Nuc__clearSpeciesCompareFunction() Description: Restore the default function to be applied during a species sort. Syntax: void Libnucnet__Nuc__clearSpeciesCompareFunction( Libnucnet__Nuc *self );Input:
Upon successful return, the species compare function for the collection has been restored to the default function. If any input is invalid, error handling is invoked. Example: Clear the species compare function in p_my_nuclei:Libnucnet__Nuc__clearSpeciesCompareFunction( p_my_nuclei ); |
Name: Libnucnet__Nuc__computeSpeciesBindingEnergy() Description: Compute the binding energy of a nuclear species in MeV. Syntax: double Libnucnet__Nuc__computeSpeciesBindingEnergy( const Libnucnet__Nuc *self, const Libnucnet__Species *p_species );Input:
Routine returns a double representing the nuclear binding energy of the input species in MeV. If input is invalid, error handling is invoked. Example: Compute the binding eneryg of oxygen-16 in the collection of nuclear species p_my_nuclei:if( ( p_species = Libnucnet__Nuc__getSpeciesByName( p_my_nuclei, "o16" ) ) ) fprintf( stdout, "Binding energy of %s = %f MeV.\n", Libnucnet__Species__getName( p_species ), Libnucnet__Nuc__computeSpeciesBindingEnergy( p_my_nuclei, p_species ) ); |
Name: Libnucnet__Nuc__extractSubset() Description: Extract a subset collection of species from another collection by an XPath expression. Syntax: Libnucnet__Nuc * Libnucnet__Nuc__extractSubset( const Libnucnet__Nuc *self const char *s_xpath );Input:
For a valid Libnucnet__Nuc structure, the routine returns a new Libnucnet__Nuc structure containing nuclei selected from the parent species collection by the XPath expression. The user must free the structure with Libnucnet__Nuc__free() when finished with it. If the input structure is not valid, Libnucnet__Nuc error handling is invoked. Examples: Extract a complete copy of Libnucnet_Nuc *p_parent and call it p_copy:p_copy = Libnucnet__Nuc__extractSubset( p_parent, NULL );Extract a subset of Libnucnet_Nuc *p_parent that includes only nuclei with Z > 10 and A < 20. p_copy = Libnucnet__Nuc__extractSubset( p_parent, "[z > 10 and a < 20]" ); |
Name: Libnucnet__Nuc__free() Description: Free the memory allocated for a Libnucnet nuclear species collection. Syntax: void Libnucnet__Nuc__free( Libnucnet__Nuc *self );Input:
Upon successful return, the memory for the collection has been freed. Example: Free up memory used for Libnucnet__Nuc *p_my_nuclei:Libnucnet__Nuc__free( p_my_nuclei ); |
Name: Libnucnet__Nuc__getLargestNucleonNumber() Description: Retrieve the largest Z, N, or A for species in the nuclide collection. Syntax: unsigned int Libnucnet__Nuc__getLargestNucleonNumber( Libnucnet__Nuc * self, const char * s_nucleon_type );Input:
Routine returns an unsigned int giving the largest Z, N, or A of the species present in the nuclide collection. If the input zone or input string for the nucleon number type is invalid, error handling is invoked. Example: Print the largest neutron number of the species present in the collection p_my_nuclei:fprintf( stdout, "The largest N is: %d\n", Libnucnet__Nuc__getLargestNucleonNumber( p_my_nuclei, "n" ) ); |
Name: Libnucnet__Nuc__getNumberOfSpecies() Description: Retrieves the number of species in the collection of nuclear species. Syntax: size_t Libnucnet__Nuc__getNumberOfSpecies( const Libnucnet__Nuc *self );Input:
Routine returns an integer representing the number of species in the collection of nuclear species. If the input is invalid, error handling is invoked. Example: Get number of species in nuclear data collection p_my_nuclei:size_t i_species_count; i_species_count = Libnucnet__Nuc__getNumberOfSpecies( p_my_nuclei ); |
Name: Libnucnet__Nuc__getSpeciesByName() Description: Retrieves the specified species. Syntax: Libnucnet__Species * Libnucnet__Nuc__getSpeciesByName( const Libnucnet__Nuc *self, const char *s_nucname );Input:
Routine returns a pointer to a Libnucnet__Species structure. If the species is not found, routine returns NULL. If the input structure is invalid, error handling is invoked. Example: Retrieve he4 from nuclear data collection p_my_nuclei:p_he4 = Libnucnet__Nuc__getSpeciesByName( p_my_nuclei, "he4" ); |
Name: Libnucnet__Nuc__getSpeciesByZA() Description: Retrieves a species in the collection of species given the Z, A, and state of the species. Syntax: Libnucnet__Species * Libnucnet__Nuc__getSpeciesByZA( const Libnucnet__Nuc *self, unsigned int i_z, unsigned int i_a, const char *s_state );Input:
Routine returns a pointer to the species retrieved. If no such species is found, the routine returns a NULL. If the input structure is not valid, error handling is invoked. Examples: Get the pointer to the species al27 in nuclear data structure p_my_nuclei:Libnucnet__Species *p_al27; p_al27 = Libnucnet__Nuc__getSpeciesByZA( p_my_nuclei, 13, 27, NULL ) );Get the pointer to species al26m in a collection of species p_my_nuclei: Libnucnet__Species *p_al26m; p_al26_m = Libnucnet__Nuc__getSpeciesByZA( p_my_nuclei, 13, 26, 'm' ) ); |
Name: Libnucnet__Nuc__getSpeciesCompareFunction() Description: Retrieve the currently set function to be applied during a species sort. Syntax: Libnucnet__Species__compare_function Libnucnet__Nuc__getSpeciesCompareFunction( Libnucnet__Nuc *self );Input:
Routine returns the currently set data the data compare function for species in the nuclide collection. If any input is invalid, error handling is invoked. Example: Set the species compare function in p_my_new_nuclei to that in p_my_old_nuclei:Libnucnet__Nuc__setSpeciesCompareFunction( p_my_new_nuclei, Libnucnet__Nuc__getSpeciesCompareFunction( p_my_old_nuclei ) ); |
Name: Libnucnet__Nuc__is_valid_input_xml() Description: Validate an input xml file for Libnucnet__Nuc. Syntax: int Libnucnet__Nuc__is_valid_input_xml( const char *s_xml_filename );Input:
For a valid input nuclear data xml file, the routine returns 1 (true). If the input file is not valid, routine returns 0 (false). For an invalid schema file, or if the schema file cannot be read over the web, routine stops and prints error message. Example: Validate the input xml file "nuclear_data.xml":if( Libnucnet__Nuc__is_valid_input_xml( "nuclear_data.xml" ) ) { printf( "Valid xml!\n" ); } |
Name: Libnucnet__Nuc__iterateSpecies() Description: Iterate through the species and apply the user-supplied iterate function. Syntax: void Libnucnet__Nuc__iterateSpecies( const Libnucnet__Nuc *self, Libnucnet__Species__iterateFunction pfFunc, void *p_user_data );Input:
The routine iterates through the species collection and applies the user-supplied routine to each species. If any input is invalid, error handling is invoked. Example: Iterate through the species in p_my_nuclei and apply the function my_iterate_function and the extra data in p_user_data:Libnucnet__Nuc__iterateSpecies( p_my_nuclei, (Libnucnet__Species__iterateFunction) my_iterate_function, p_user_data ); |
Name: Libnucnet__Nuc__new() Description: Create a new Libnucnet__Nuc structure. Syntax: Libnucnet__Nuc *Libnucnet__Nuc__new( );Output: The routine returns a pointer to a new Libnucnet__Nuc structure, that is, a collection of nuclear species. If it is not possible to allocate memory for the new structure, Libnucnet error handling is invoked. Example: Create a Libnucnet__Nuc structure p_my_nuclei:p_my_nuclei = Libnucnet__Nuc__new( ); |
Name: Libnucnet__Nuc__new_from_xml() Description: Creates a Libnucnet__Nuc structure from an input xml nuclear data file. Syntax: Libnucnet__Nuc * Libnucnet__Nuc__new_from_xml( const char *s_xml_filename, const char *s_xpath );Input:
For a valid input nuclear data xml file, the routine returns a pointer to a new Libnucnet__Nuc structure. If the routine cannot allocate enough memory or if the xpath expression is invalid, Libnucnet__Nuc error handling is invoked. Examples: Store the nuclear data in "nuclear_data.xml" to p_my_nuclei:p_my_nuclei = Libnucnet__Nuc__new_from_xml( "nuclear_data.xml", NULL );Store the nuclear data for neutrons and protons and nuclei with z >= 6 in "nuclear_data.xml" to p_my_nuclei: p_my_nuclei = Libnucnet__Nuc__new_from_xml( "nuclear_data.xml", "[ a = 1 or z >= 6 ]" ); |
Name: Libnucnet__Nuc__new_from_xml_text_reader() Description: Reads in nuclear network nuclear data from xml file and stores the species in a structure. This routine does not build a tree, so the memory overhead is greatly reduced. The species to be stored are selected with an xpath expression, but unlike the routine that builds with a tree, the xpath expression is relative to the individual species. Syntax: Libnucnet__Nuc * Libnucnet__Nuc__new_from_xml( const char *s_xml_filename, const char *s_xpath );Input:
If the input is valid, routine returns a pointer to a Libnucnet species structure containing species selected by the xpath expression. If no species were found, p_species_list is empty. If the input file or xpath expression(s) is invalid, Libnucnet__Nuc error handling is invoked. Examples: Store all the species from "nuclides.xml" to p_species_list:p_species_list = Libnucnet__Nuc__new_from_xml( "nuclides.xml", NULL );Store the species from "nuclides.xml" with z < 20 and a < 100 in p_species_list: p_species_list = Libnucnet__Nuc__new_from_xml( "nuclides.xml", ".//z[. < 20 and (preceding-sibling::a[. < 100] or following-sibling::a[. < 100])]" ); |
Name: Libnucnet__Nuc__removeSpecies() Description: Remove a species from a Libnucnet__Nuc structure. Syntax: int Libnucnet__Nuc__removeSpecies( Libnucnet__Nuc *self, Libnucnet__Species *p_species );Input:
Routine returns 1 for successful removal, 0 for failure. Upon successful return, the species has been removed from the Libnucnet__Nuc structure and the species indices have been updated. If the input is invalid, Libnucnet__Nuc error handling is invoked. Example: Remove he4 from Libnucnet__Nuc *p_my_nuclei:Libnucnet__Species *p_he4 = Libnucnet__Nuc__getSpeciesByName( p_my_nuclei, "he4" ); if( Libnucnet__Nuc__removeSpecies( p_my_nuclei, p_he4 ) ) printf( "Successful removal.\n" ); |
Name: Libnucnet__Nuc__setSpeciesCompareFunction() Description: Set the function to be applied during a species sort. Syntax: void Libnucnet__Nuc__setSpeciesCompareFunction( Libnucnet__Nuc *self, Libnucnet__Species__compare_function pfFunc );Input:
Upon successful return, the data compare function for species in the collection has been set to the input function. If any input is invalid, error handling is invoked. Example: Set the species compare function in p_my_nuclei to my_sort_function:Libnucnet__Nuc__setSpeciesCompareFunction( p_my_nuclei, (Libnucnet__Species__compare_function) my_sort_function ); |
Name: Libnucnet__Nuc__sortSpecies() Description: Sort the species in a Libnucnet__Nuc structure according to the current Libnucnet__Species__compare_function. Syntax: void Libnucnet__Nuc__sortSpecies( Libnucnet__Nuc *self );Input:
On successful return, the species in the structure have been sorted according to the current Libnucnet__Species__compare_function. If the default is set, the species will be sorted in increasing Z and A. If the input structure is not valid or memory cannot be allocated for the sorting, error handling is invoked. Example: Sort the species in p_my_nuclei:Libnucnet__Nuc__sortSpecies( p_my_nuclei ); |
Name: Libnucnet__Nuc__updateFromXml() Description: Updates the data in a Libnucnet__Nuc structure from an input xml nuclear data file. Syntax: void Libnucnet__Nuc__updateFromXml( Libnucnet__Nuc *self, const char *s_xml_filename, const char *s_xpath );Input:
For a valid Libnucnet__Nuc structure and a valid input nuclear data xml file, the routine updates the structure with the data in the input data file. If the routine cannot allocate enough memory or if the xpath expression is invalid, Libnucnet__Nuc error handling is invoked. Example: Update the data in p_my_nuclei with data in "new_nuclear_data.xml":Libnucnet__Nuc__updateFromXml( p_my_nuclei, "new_nuclear_data.xml", NULL ); |
Name: Libnucnet__Nuc__updateFromXmlTextReader() Description: Updates nuclear network nuclide data in a Libnucnet__Nuc structure from and xml file. This routine does not build a tree, so the memory overhead is greatly reduced. The species to be stored are selected with an xpath expression, but unlike the routine that builds with a tree, the xpath expression is relative to the individual species. Syntax: void Libnucnet__Nuc__updateFromXmlTextReader( Libnucnet__Nuc *self, const char *s_xml_filename, const char *s_xpath );Input:
If the input is valid, the input Libnucnet__Nuc structure is updated with the species data from the input file. Example: Update the species data in p_my_species, with data from the file "nuclides_new.xml"Libnucnet__Nuc__updateFromXmlTextReader( p_my_species, "nuclides_new.xml", NULL ); |
Name: Libnucnet__Nuc__updateSpecies() Description: Update a species in a species collection. Syntax: int Libnucnet__Nuc__updateSpecies( Libnucnet__Nuc *self, Libnucnet__Species *p_species );Input:
Routine returns 1 for successful update and 0 for failure. For valid input, the routine updates the Libnucnet__Nuc structure with the new species. If the species already exists, its data are replaced with the new data. If not, the species is added. Example: Update p_species in the species collection p_my_nuclei:if( Libnucnet__Nuc__updateSpecies( p_my_nuclei, p_species ) ) printf( "Species update successful\n." ); |
Name: Libnucnet__Nuc__writeToXmlFile() Description: Output a Libnucnet__Nuc structure to an xml file. Syntax: void Libnucnet__Nuc__writeToXmlFile( const Libnucnet__Nuc *self, const char *s_xml_filename );Input:
Upon successful return, the contents of the collection of species have been written to an xml file. If the input is invalid, Libnucnet__Nuc error handling is invoked. Example: Dump the contents of Libnucnet__Nuc *p_my_nuclei to the xml file my.xml:Libnucnet__Nuc__writeToXmlFile( p_my_nuclei, "my.xml" ); |
Name: Libnucnet__Nuc__writeToXmlFileWithTextWriter() Description: Output a Libnucnet__Nuc structure to an xml file. Use the text writer which can save on memory usage. Syntax: void Libnucnet__Nuc__writeToXmlFileWithTextWriter( const Libnucnet__Nuc *self, const char *s_xml_filename int i_indent );Input:
Upon successful return, the contents of the collection of species have been written to an xml file. If the input is invalid, Libnucnet__Nuc error handling is invoked. Example: Dump the contents of Libnucnet__Nuc *p_my_nuclei to the xml file my.xml and indent the output:Libnucnet__Nuc__writeToXmlFileWithTextWriter( p_my_nuclei, "my.xml", 1 ); |
Name: Libnucnet__Species__computePartitionFunction() Description: Compute the partition function of a species in the collection of nuclear species given the temperature. Syntax: double Libnucnet__Species__computePartitionFunction( const Libnucnet__Species *self, double d_t9 );Input:
Routine returns a double representing the partition function of the species retrieved. If input is invalid, error handling is invoked. Example: Compute the partition function of he4 in the collection of nuclear species p_my_nuclei at a t9 of 0.2:if( ( p_species = Libnucnet__Nuc__getSpeciesByName( p_my_nuclei, "he4" ) ) ) { d_partf = Libnucnet__Species__computePartitionFunction( p_species, 0.2 ); } |
Name: Libnucnet__Species__computeQuantumAbundance() Description: Compute the quantum abundance of a nuclear species at the given the temperature and density Syntax: double Libnucnet__Species__computeQuantumAbundance( const Libnucnet__Species *self, double d_t9, double d_rho );Input:
Routine returns a double representing the quantum abundance of the given species. If input is invalid, error handling is invoked. Example: Compute the quantum abundance of si29 in the collection of nuclear species p_my_nuclei at a t9 of 0.2 and mass density of 1000 g/cc:if( ( p_species = Libnucnet__Nuc__getSpeciesByName( p_my_nuclei, "si29" ) ) ) { d_y_Q = Libnucnet__Species__computeQuantumAbundance( p_species, 0.2, 1000. ); } |
Name: Libnucnet__Species__copy() Description: Copy a Libnucnet__Species. Syntax: Libnucnet__Species * Libnucnet__Species__copy( const Libnucnet__Species *p_species );Input:
The routine returns a pointer to a new Libnucnet__Species structure that is a copy of the input structure. If it is not possible to allocate memory for the new structure, Libnucnet error handling is invoked. Example: Create a copy of p_species:p_copy = Libnucnet__Species__copy( p_species ); |
Name: Libnucnet__Species__createLatexString() Description: Create a string giving the name of the species in latex format. Syntax: char * Libnucnet__Species__createLatexString( const Libnucnet__Species * self );Input:
Routine returns a new string giving the name of the species in latex form. If the input species is invalid, error handling is invoked. Example: Output the name of species p_species in latex form:s_latex_name = Libnucnet__Species__createLatexString( p_species ); fprintf( stdout, "Here is the latex form of the species name (with math delimiters): $%s$.\n" s_latex_name ); free( s_latex_name ); |
Name: Libnucnet__Species__free() Description: Free the memory allocated for a Libnucnet species. Syntax: void Libnucnet__Species__free( Libnucnet__Species *self );Input:
Upon successful return, the memory for the species has been freed. Example: Free up memory used for Libnucnet__Species *p_species:Libnucnet__Species__free( p_species ); |
Name: Libnucnet__Species__getA() Description: Retrieves the A (that is, the mass number) of a species in the collection of nuclear species. Syntax: unsigned int Libnucnet__Species__getA( const Libnucnet__Species *self );Input:
Routine returns the A of the input species. If the species is invalid, Libnucnet__Nuc error handling is invoked. Example: Print the A of au197 in p_my_nuclei (should of course be 197!):p_species = Libnucnet__Nuc__getSpeciesByName( p_my_nuclei, "au197" ); if( p_species ) printf( "A = %d\n", Libnucnet__Species__getA( p_species ) ); |
Name: Libnucnet__Species__getIndex() Description: Routine to return the index of a species in a Libnucnet__Nuc structure. Syntax: size_t Libnucnet__Species__getIndex( const Libnucnet__Species *self );Input:
Returns the index of the species in the nuclear species collection. If the species is invalid, Libnucnet__Nuc error handling is invoked. Example: Print the index of he4 in the Libnucnet structure p_my_nuclei:if( p_he4 = Libnucnet__Nuc__getSpeciesByName( p_my_nuclei, "he4" ) ) { printf( "Index of he4 = %d\n", Libnucnet__Species__getIndex( p_he4 ) ); } |
Name: Libnucnet__Species__getMassExcess() Description: Retrieves the mass excess of a species in the nuclear species collection. Syntax: double Libnucnet__Species__getMassExcess( const Libnucnet__Species *self );Input:
Routine returns a double giving the mass excess of the requested species. If the input species is invalid, error handling is invoked. Example: Print the mass excess of sn120 in the structure Libnucnet__Nuc *p_my_nuclei:p_species = Libnucnet__Nuc__getSpeciesByName( p_my_nuclei, "sn120" ); if( p_species ) { printf( "Mass excess (MeV) = %lf\n, Libnucnet__Species__getMassExcess( p_species ) ); } |
Name: Libnucnet__Species__getName() Description: Retrieves the name of a species in the collection of nuclear species. Syntax: const char * Libnucnet__Species__getName( const Libnucnet__Species *self );Input:
Routine returns a string representing the name of the requested species. If the species is not valid, Libnucnet__Nuc error handling is invoked. Example: Print the name of the Libnucnet__Species *p_species:printf( "Name of species = %s\n", Libnucnet__Species__getName( p_species ) ); |
Name: Libnucnet__Species__getPartitionFunctionLog10() Description: Retrieves the pointer to the gsl_vector containing the data for the partition function for a species in the nuclear species collection. Syntax: gsl_vector * Libnucnet__Species__getPartitionFunctionLog10( const Libnucnet__Species *self );Input:
Routine returns the pointer to the gsl_vector containing the log10 partition function data points for the species. If the input species is invalid, error handling is invoked. Example: Print the log10 partition function data for mo94 in *p_my_nuclei:p_species = Libnucnet__Nuc__getSpeciesByName( p_my_nuclei, "mo94" ); if( p_species ) gsl_vector_fprintf( stdout, Libnucnet__Species__getPartitionFunctionLog10( p_species ), "%g" ); |
Name: Libnucnet__Species__getPartitionFunctionT9() Description: Retrieves the pointer to the gsl_vector containing the temperature data for the partition function for a species in the nuclear species collection. Syntax: gsl_vector * Libnucnet__Species__getPartitionFunctionT9( const Libnucnet__Species *self );Input:
Routine returns the pointer to the gsl_vector containing the T9 points for the partition function data for the species. If the input species is invalid, error handling is invoked. Example: Print the gsl_vector t9 array for mo94 in *p_my_nuclei:p_species = Libnucnet__Nuc__getSpeciesByName( p_my_nuclei, "mo94" ); if( p_species ) gsl_vector_fprintf( stdout, Libnucnet__Species__getPartitionFunctionT9( p_species ), "%g" ); |
Name: Libnucnet__Species__getSource() Description: Retrieves the message about the source of data for a species in the collection of nuclear species. Syntax: const char * Libnucnet__Species__getSource( const Libnucnet__Species *self );Input:
Routine returns a string containing information about the data for the requested species. If there are no data, an empty string is returned. If the species is not valid, Libnucnet__Nuc error handling is invoked. Example: Print the source information about the Libnucnet__Species *p_species:printf( "Data source for species = %s\n", Libnucnet__Species__getSource( p_species ) ); |
Name: Libnucnet__Species__getSpin() Description: Retrieves the spin of a species in the collection of nuclear species. Syntax: double Libnucnet__Species__getSpin( const Libnucnet__Species *self );Input:
Routine returns a double giving the spin of the requested species. If the species is invalid, Libnucnet__Nuc error handling is invoked. Example: Print the spin of bi209 in the Libnucnet__Nuc structure *p_my_nuclei:if( ( p_species = Libnucnet__Nuc__getSpeciesByName( p_my_nuclei, "bi209" ) ) ) { printf( "Spin = "%lf\n", Libnucnet__Species__getSpin( p_species ) ); } |
Name: Libnucnet__Species__getStateString() Description: Retrieves the state for a species in the collection of nuclear species. Syntax: const char * Libnucnet__Species__getStateString( const Libnucnet__Species *self );Input:
Routine returns a string containing the state for the requested species. If the species has no state string (for example, the species has only one state), the routine returns NULL. If the species is not valid, Libnucnet__Nuc error handling is invoked. Example: Print the state of the Libnucnet__Species *p_species:if( Libnucnet__Species__getStateString( p_species ) ) { printf( "State for species = %s\n", Libnucnet__Species__getStateString( p_species ) ); } |
Name: Libnucnet__Species__getZ() Description: Retrieves the Z (that is, the atomic number) of a species in the collection of nuclear species. Syntax: unsigned int Libnucnet__Species__getZ( const Libnucnet__Species *self );Input:
Routine returns the Z of the input species. If the species is invalid, Libnucnet__Nuc error handling is invoked. Example: Print the Z of au197 in p_my_nuclei:p_species = Libnucnet__Nuc__getSpeciesByName( p_my_nuclei, "au197" ); if( p_species ) printf( "Z = %d\n", Libnucnet__Species__getZ( p_species ) ); |
Name: Libnucnet__Species__new() Description: Create a new species from the input data. Syntax: Libnucnet__Species * Libnucnet__Species__new( unsigned int i_z, unsigned int i_a, const char *s_source, int i_state_flag, const char *s_state, double d_mass_excess, double d_spin, gsl_vector *p_t9, gsl_vector *p_log10_partf );Input:
For valid input, the routine returns a pointer to a new Libnucnet__Species structure. If the species cannot be created, error handling is invoked. Examples: Create li7 (Z = 3, A = 7, spin = 3/2, mass excess = 14.908 MeV), and partition function data stored in gsl_vectors p_t9 and p_log10_partf:p_li7 = Libnucnet__Species__new( 3, 7, "example", 0, "", 14.908, 1.5, p_t9, p_log10_partf );Create ground state (g) of al26 (Z = 13, A = 26, spin = 5, mass excess = -12.210 MeV), and partition function data stored in p_t9_g and p_log10_partf_g and metastable state (m) of al26 (Z = 13, A = 26, spin = 0, mass excess = -11.982), and partition function data stored in p_t9_m and p_log10_partf_m: p_al26g = Libnucnet__Species__new( 13, 26, "example data", 1, "g", -12.210, 5., p_t9_g, p_log10_partf_g ); p_al26m = Libnucnet__Species__new( 13, 26, "example data", 1, "m", -11.982, 0., p_t9_m, p_log10_partf_m ); |
Name: Libnucnet__Species__updateMassExcess() Description: Update the mass excess of the species. Syntax: void Libnucnet__Species__updateMassExcess( Libnucnet__Species *self, double d_new_mass_excess );Input:
Upon successful return, the mass excess for the species has been updated to the new value. If the input species is invalid, error handling is invoked. Example: Update the mass excess of p_species to 8.22 MeV:Libnucnet__Species__updateMassExcess( p_species, 8.22 ); |
Name: Libnucnet__Species__updatePartitionFunctionData() Description: Update the partition function data for the species. Syntax: void Libnucnet__Species__updatePartitionFunctionData( Libnucnet__Species *self, gsl_vector *p_new_t9_array, gsl_vector *p_new_log10_partf_array );Input:
Upon successful return, the partition function data for the species has been updated to the new values. If the input species is invalid, or if the new arrays cannot be allocated, error handling is invoked. Example: Update the t9 and log10 partf arrays for the species p_species to p_new_t9 and p_new_log10_partf:Libnucnet__Species__updatePartitionFunctionData( p_species, p_new_t9, p_new_log10_partf ); |
Name: Libnucnet__Species__updateSource() Description: Update the string giving the data source for the species. Syntax: void Libnucnet__Species__updateSource( Libnucnet__Species *self, const char *s_source );Input:
Upon successful return, the data source string for the species has been updated to the new value. If the input species is invalid, error handling is invoked. Example: Update the source string for the data to "Updated data source":Libnucnet__Species__updateSource( p_species, "Updated data source" ); |
Name: Libnucnet__Species__updateSpin() Description: Update the spin of the species. Syntax: void Libnucnet__Species__updateSpin( Libnucnet__Species *self, double d_new_spin );Input:
Upon successful return, the spin of the species has been updated to the new value. If the input species is invalid, error handling is invoked. Example: Update the spin of p_species to 3/2:Libnucnet__Species__updateSpin( p_species, 1.5 ); |