Structures |
---|
Name: Libnucnet__Net Description: Libnucnet__Net is a structure that stores data about the nuclei and reactions among them in a nuclear reaction network. A Libnucnet__Net structure comprises a Libnucnet__Nuc structure and a Libnucnet__Reac structure. The contents of Libnucnet__Net are not made public by the API but rather are accessed through the API routines. |
Name: Libnucnet__NetView Description: Libnucnet__NetView is a structure that stores a view of Libnucnet__Net structure. A view is a subset of the parent Libnucnet__Net. The nuclei included are chosen by XPath expressions as are the reactions, with the constraint that the reactions must be valid for the given selection of nuclei. The view structure contains one element, a Libnucnet__Net, which may be passed into all API routines that take such structures as input. It is important to note that the Libnucnet__Net member of a view does not own the nuclide or reaction data, so the user must be careful in modifying those data. |
Name: Libnucnet__Zone Description: Libnucnet__Zone is a structure that stores data about the abundances of nuclear species and the reaction rates among them in a particular zone. A zone is labelled by three strings. If only one zone exists, the default labels are "0", "0", "0"; however, the user may choose other labels for the single zone. A Libnucnet__Zone structure also contains a pointer to the Libnucnet__Net structure. A Libnucnet__Zone structure thus contains data that can change with each timestep in a calculation (abundances and reaction rates) along with data that is essentially fixed (Libnucnet__Net data). The contents of a Libnucnet__Zone structure are not made public by the API. |
Name: Libnucnet Description: Libnucnet is a structure that stores data about nuclear reaction networks. It is in fact composed of a Libnucnet__Net structure and a hash of Libnucnet__Zone structures. The contents of a Libnucnet structure are not made public by the API but rather are accessed through the API routines. |
User-Supplied Routines |
---|
Name: Libnucnet__NetView__iterateFunction() Description: User-supplied routine to be applied during an iteration over the network views in a zone. Syntax: void Libnucnet__NetView__iterateFunction( Libnucnet__NetView * p_view, const char *s_label1, const char *s_label2, const char *s_label3, void *p_user_data );Input:
On successful return, the function has been applied. |
Name: Libnucnet__Zone__compare_function() Description: User-supplied routine to sort the zones during a zone iteration. Syntax: int Libnucnet__Zone__compare_function( const Libnucnet__Zone *p_zone1, const Libnucnet__Zone *p_zone2 );Input:
User's routine must return -1 if zone 1 is less than zone 2, 0 if the two zones are equal, and 1 if zone 1 is greater than zone 2. If this routine is not supplied through Libnucnet__setZoneCompareFunction, the default is not to sort the zones but rather to iterate them in the order in which they are stored internally, which is not defined by the user. |
Name: Libnucnet__Zone__iterateFunction() Description: User-supplied routine to be applied during an iteration over the zones. Syntax: int Libnucnet__Zone__iterateFunction( Libnucnet__Zone *self, void *p_user_data );Input:
User's routine must return 1 to continue or 0 to stop. |
Name: Libnucnet__Zone__optional_property_iterate_function() Description: User-supplied routine to be applied during an iteration over the optional properties in a zone. Syntax: void Libnucnet__Zone__optional_property_iterate_function( const char *s_name, const char *s_tag1, const char *s_tag2, const char *s_value, void *p_user_data );Input:
On successful return, the function has been applied. |
Name: Libnucnet__Zone__screeningFunction() Description: Optional user-supplied screening function. Syntax: void Libnucnet__Zone__screeningFunction( Libnucnet__Zone * self, Libnucnet__Reaction * p_reaction, double d_t9, double d_rho, double d_ye, double * p_forward_rate, double * p_reverse_rate );Input:
User's routine must compute the screening for a reaction at the input temperature, density, and Ye and apply the screening to the forward and reverse rates for the reaction. The user can set this routine with Libnucnet__Zone__setScreeningFunction for libnucnet to use from Libnucnet__Zone__computeRates(). |
Routines |
---|
Name: Libnucnet__NetView__addReaction() Description: Add a reaction to a network view. Syntax: int Libnucnet__NetView__addReaction( Libnucnet__NetView * self, Libnucnet__Reaction * p_reaction );Input:
Returns 1 (true) if the addition succeeded or 0 (false) if not (either the reaction was not valid for the view or could not be added). If input is invalid, Libnucnet error handling is invoked. Example: Add reaction p_reaction to network view p_view:if( Libnucnet__NetView__addReaction( p_view, p_reaction ) ) fprintf( stdout, "Addition succeeded.\n" ); |
Name: Libnucnet__NetView__copy() Description: Get a copy of a network view. Syntax: Libnucnet__NetView * Libnucnet__NetView__copy( Libnucnet__NetView * self );Input:
Returns a new network view that is a copy of the input network view. If the input structure is invalid, Libnucnet error handling is invoked. Example: Copy p_view:p_new_view = Libnucnet__NetView__copy( p_view ); |
Name: Libnucnet__NetView__free() Description: Free the memory for a Libnucnet__NetView structure. Syntax: void Libnucnet__NetView__free( Libnucnet__NetView *self );Input:
On return, the memory allocated for the network view has been cleared. If the input pointer is not valid, Libnucnet error handling is invoked. Example: Free the memory for the network view p_my_view:Libnucnet__NetView__free( p_my_view ); |
Name: Libnucnet__NetView__getNet() Description: Retrieve the Libnucnet__Net in a Libnucnet__NetView structure. Syntax: Libnucnet__Net * Libnucnet__NetView__getNet( Libnucnet__NetView *self );Input:
Routine returns the pointer to the Libnucnet__Net structure in a valid network view. If the input pointer is not valid, Libnucnet error handling is invoked. Example: Get the network in view p_my_view:p_net = Libnucnet__NetView__getNet( p_my_view ); |
Name: Libnucnet__NetView__new() Description: Create a new Libnucnet__NetView structure. Syntax: Libnucnet__NetView * Libnucnet__NetView__new( Libnucnet__Net * p_net, const char * s_nuc_xpath, const char * s_reac_xpath );Input:
Routine returns a new network view. It is the user's responsibility to free the view once it is no longer needed. If the input is not valid, Libnucnet error handling is invoked. Example: Create a network view of the parent network p_net containing all nuclei and all valid reactions among them:p_view = Libnucnet__NetView( p_net, " ", " " ); |
Name: Libnucnet__NetView__removeReaction() Description: Remove a reaction from a network view. Syntax: int Libnucnet__NetView__removeReaction( Libnucnet__NetView * self, Libnucnet__Reaction * p_reaction );Input:
Returns 1 (true) if the removal succeeded or 0 (false) if not. If input is invalid, Libnucnet error handling is invoked. Example: Remove reaction p_reaction from network view p_view:if( Libnucnet__NetView__removeReaction( p_view, p_reaction ) ) fprintf( stdout, "Removal succeeded.\n" ); |
Name: Libnucnet__NetView__wasNetUpdated() Description: Check if the underlying network of a network view has been updated since the view was generated. Syntax: int Libnucnet__NetView__wasNetUpdated( Libnucnet__NetView * self );Input:
Returns 1 (true) if the underlying network has been updated or 0 (false) if not. If the network has been updated, the user should generate a new view. If input is invalid, Libnucnet error handling is invoked. Example: Check whether the underlying network of p_view has been updated:if( Libnucnet__NetView__wasNetUpdated( p_view ) ) fprintf( stdout, "The network has been updated.\n" ); |
Name: Libnucnet__Net__computeRatesForReaction() Description: Compute the forward and reverse rates for a valid reaction at the input temperature and density. The rates are those between reacting species. Syntax: void Libnucnet__Net__computeRatesForReaction( Libnucnet__Net *self, Libnucnet__Reaction *p_reaction, double d_t9, double d_rho, void *p_extra_data, double *p_forward, double *p_reverse );Input:
On successful return, the routine returns the computed forward and reverse rates. If any input is invalid, error handling is invoked. Example: Compute the forward and reverse rates for c12 + h1 -> n13 + gamma at a temperature T9 = 1 and density rho = 2.e6 g/cc in Libnucnet__Net *p_my_nucnet:p_reaction = Libnucnet__Reac__getReactionByString( Libnucnet__Net__getReac( p_my_nucnet ), "c12 + h1 -> n13 + gamma" ); Libnucnet__Net__computeRatesForReaction( p_my_nucnet, p_reaction, 1., 2.e6, NULL, &d_forward, &d_reverse ); |
Name: Libnucnet__Net__computeReactionQValue() Description: Compute the Q value (in MeV) for the specified reaction in the Libnucnet__Net structure. Syntax: double Libnucnet__Net__computeReactionQValue( Libnucnet__Net *self, Libnucnet__Reaction *p_reaction, );Input:
Routine returns a double containing the Q value of the reaction. If either the Net or Reaction input structure is invalid, error handling is invoked. Example: Compute the Q value of the reaction ne20 + he4 -> mg24 + gamma from Libnucnet__Net *p_my_net and store the value in d_Qvalue:d_Qvalue = Libnucnet__Net__computeReactionQValue( p_my_net, Libnucnet__Reac__getReactionByString( Libnucnet__Net__getReac( p_my_net ), "ne20 + he4 -> mg24 + gamma" ) ); |
Name: Libnucnet__Net__computeReverseRatioCorrectionFactorForReaction() Description: Use the user-supplied NSE correction factor function and its associated data to compute the reverse ratio correction factor for a reaction. Syntax: double Libnucnet__Net__computeReverseRatioCorrectionFactorForReaction( Libnucnet__Net *self, Libnucnet__Reaction *p_reaction, double d_t9, double d_rho, double d_ye, Libnucnet__Species__nseCorrectionFactorFunction p_user_function, void *p_user_data );Input:
For valid input, the routine returns a double with the reverse ratio correction factor as computed by the user's supplied function and data. Examples: Compute the reverse ratio correction factor for the reaction c12 + h1 -> n13 + gamma at t9 = 2, rho = 1.e11, and Ye = 0.5. Use the user-supplied screening routine my_correction_function and pass no extra data for the correction calculation:p_reaction = Libnucnet__Reac__getReactionByString( Libnucnet__Net__getReac( p_my_net ), "c12 + h1 -> n13 + gamma" ); d_correction = Libnucnet__Net__computeReverseRatioCorrectionFactorForReaction( p_my_net, p_reaction, 2., 1.e11, 0.5, (Libnucnet__Species__nseCorrectionFactorFunction) my_correction_function, NULL );Compute the reverse ratio correction factor for the reaction c12 + h1 -> n13 + gamma at t9 = 2, rho = 1.e11, and Ye = 0.5. Use the user-supplied correction routine my_correction_function and pass extra data for the correction calculation through the data structure my_data: struct my_data = { double x; double y; }; my_data.x = 1.; my_data.y = 2.; p_reaction = Libnucnet__Reac__getReactionByString( Libnucnet__Net__getReac( p_my_net ), "c12 + h1 -> n13 + gamma" ); d_correction = Libnucnet__Net__computeReverseRatioCorrectionFactorForReaction( p_my_net, p_reaction, 2., 1.e11, 0.5, (Libnucnet__Species__nseCorrectionFactorFunction) my_correction_function, &my_data ); |
Name: Libnucnet__Net__createValidReactionLatexString() Description: Create a latex string for a valid reaction. Syntax: char * Libnucnet__Net__createValidReactionLatexString( const Libnucnet__Net * self, const Libnucnet__Reaction * p_reaction );Input:
Routine returns a new string for the reaction in latex form. If the reaction is invalid, routine returns NULL. If the input network or reaction pointer is invalid, error handling is invoked. Example: Output the latex string for p_reaction in p_net:s_latex_string = Libnucnet__Net__createValidReactionLatexString( p_net, p_reaction ); fprintf( stdout, "Here is the latex form of the reaction (with math delimiters): $%s$.\n" s_latex_string ); free( s_latex_string ); |
Name: Libnucnet__Net__free() Description: Free the memory for a Libnucnet__Net structure. Syntax: void Libnucnet__Net__free( Libnucnet__Net *self );Input:
On return, the memory allocated for the network has been cleared. If the input pointer is not valid, Libnucnet error handling is invoked. Example: Free the memory for the network in p_my_nucnet:Libnucnet__Net__free( p_my_nucnet ); |
Name: Libnucnet__Net__getNuc() Description: Retrieve the Libnucnet__Nuc structure from a Libnucnet__Net structure. Syntax: Libnucnet__Nuc * Libnucnet__Net__getNuc( const Libnucnet__Net *self );Input:
Routine returns the pointer to the Libnucnet__Nuc structure contained in the Libnucnet__Net structure. If the input pointer is not valid, Libnucnet error handling is invoked. Example: Retrieve the nuclear data structure p_my_nuclei from the network structure p_my_network:p_my_nuclei = Libnucnet__Net__getNuc( p_my_network ); |
Name: Libnucnet__Net__getNumberOfValidReactions() Description: Routine to count the number of valid reactions, that is, those reactions that are between species within the network. Exclude those that are not from use in the network. Also check reactions for conservation of nucleon number, charge, and lepton number. Syntax: size_t Libnucnet__Net__getNumberOfValidReactions( Libnucnet__Net *self );Input:
Routine returns the number of valid reactions. If the input network structure is not valid, Libnucnet__Net error handling is invoked. Example: Print the number of valid reactions in network structure p_my_network:printf( "Number of valid reactions = %lu\n", (unsigned long) Libnucnet__Net__getNumberOfValidReactions( p_my_network ) ); |
Name: Libnucnet__Net__getReac() Description: Retrieve the Libnucnet__Reac structure from a Libnucnet__Net structure. Syntax: Libnucnet__Reac * Libnucnet__Net__getReac( const Libnucnet__Net *self );Input:
Routine returns the pointer to the Libnucnet__Reac structure contained in the Libnucnet__Net structure. If the input pointer is not valid, Libnucnet error handling is invoked. Example: Retrieve the nuclear reaction data structure p_my_reactions from the network structure p_my_network:p_my_reactions = Libnucnet__Net__getReac( p_my_network ); |
Name: Libnucnet__Net__isValidReaction() Description: Determine whether a reaction is valid, that is, that it is between species in the network and conserves charge and baryon and lepton number. Syntax: int Libnucnet__Net__isValidReaction( const Libnucnet__Net *self, const Libnucnet__Reaction *p_reaction );Input:
For valid input, routine returns 1 if the reaction is valid and 0 if it is not. Example: Determine whether the reaction "o18 + he4 -> ne22 + gamma" is valid in the network Libnucnet__Net *p_my_net, and if so, print it out:p_reaction = Libnucnet__Reac__getReactionByString( Libnucnet__Net__getReac( p_my_net ), "o18 + he4 -> ne22 + gamma" ); if( Libnucnet__Net__isValidReaction( p_my_net, p_reaction ) ) { printf( "%s is a valid reaction!\n", Libnucnet__Reaction__getString( p_reaction ) ); } |
Name: Libnucnet__Net__is_valid_input_xml() Description: Validate an xml file for Libnucnet__Net. Syntax: int Libnucnet__Net__is_valid_input_xml( const char *s_xml_filename );Input:
For a valid input nuclear network data xml file, the routine returns 1 (true). For an invalid file, routine returns 0 (false). If the schema file is invalid, or if the schema file cannot be read over the web, routine stops and prints error message. Example: Validate the input Libnucnet__Net xml file "network_data.xml":if( Libnucnet__Net__is_valid_input_xml( "network_data.xml" ) ) printf( "Valid xml!\n" ); |
Name: Libnucnet__Net__new() Description: Creates a new Libnucnet__Net structure. Syntax: Libnucnet__Net *Libnucnet__Net__new( );Output: Routine returns a pointer to a new Libnucnet__Net structure. If the new structure cannot be created, error handling is invoked. Example: Create a new network structure p_my_network:p_my_network = Libnucnet__Net__new( ); |
Name: Libnucnet__Net__new_from_xml() Description: Reads in nuclear and reaction data from an xml file and stores the reactions in a Libnucnet__Net structure. The species and the reactions to be stored may be selected with an xpath expression. Syntax: Libnucnet__Net * Libnucnet__Net__new_from_xml( const char *s_xml_filename, const char *s_xpath_nuc, const char *s_xpath_reac );Input:
If the input is valid, routine returns a pointer to a Libnucnet__Net structure containing data selected by the xpath expressions. If no data were found, the return is NULL. Examples: Store all the data from "input_data.xml" to p_my_network:p_my_network = Libnucnet__Net__new_from_xml( "input_data.xml", NULL, NULL );Store the data from "input_data.xml", but only those nuclides with Z less than 4 and only reactions containing 'h1' as a reactant, to p_my_network: p_my_network = Libnucnet__Net__new_from_xml( "input_data.xml", "[z < 4]", "[reactant = 'h1']" ); |
Name: Libnucnet__Net__updateFromXml() Description: Reads in nuclear and reaction data from xml file and updates the data in the input Libnucnet__Net structure. The nuclides and reactions to be stored may be selected with an xpath expression. Syntax: void Libnucnet__Net__updateFromXml( Libnucnet__Net *self, const char *s_xml_filename, const char *s_xpath_nuc, const char *s_xpath_reac );Input:
If the input is valid, routine updates data in the structure with that contained in the input file. If the input structure is invalid or the input file does not exist, error handling is invoked. Examples: Update the data in Libnucnet__Net *p_my_network with data from "input_new_data.xml":Libnucnet__Net__new_from_xml( p_my_network, "input_new_data.xml", NULL, NULL );Update the data in Libnucnet__Net *p_my_network with data for nuclides with Z less than 10 from "input_new_data.xml": Libnucnet__Net__updateFromXml( p_my_network, "input_new_data.xml", "[z < 10]", NULL ); |
Name: Libnucnet__Net__updateFromXmlTextReader() Description: Reads in nuclear and reaction data from xml file and updates the data in the input Libnucnet__Net structure without building a tree. The nuclides and reactions to be stored may be selected with an xpath expression. In this routine, xpath expressions are relative to individual nuclides and reactions. Syntax: void Libnucnet__Net__updateFromXmlTextReader( Libnucnet__Net *self, const char *s_xml_filename, const char *s_xpath_nuc, const char *s_xpath_reac );Input:
If the input is valid, routine updates data in the structure with that contained in the input file. If the input structure is invalid or the input file does not exist, error handling is invoked. Examples: Update the data in Libnucnet__Net *p_my_network with data from "input_new_data.xml":Libnucnet__Net__new_from_xml( p_my_network, "input_new_data.xml", NULL, NULL );Update the data in Libnucnet__Net *p_my_network with data for nuclides with Z less than 10 from "input_new_data.xml": Libnucnet__Net__updateFromXml( p_my_network, "input_new_data.xml", ".//z[. < 10]", NULL ); |
Name: Libnucnet__Net__writeToXmlFile() Description: Output a Libnucnet__Net structure to an xml file. Syntax: void Libnucnet__Net__writeToXmlFile( const Libnucnet__Net *self, const char *s_xml_filename );Input:
Upon successful return, the contents of the network are written to s_xml_filename. Example: Dump the contents of Libnucnet__Net *p_my_network to the xml file my.xml:Libnucnet__Net__writeToXmlFile( p_my_network, "my.xml" ); |
Name: Libnucnet__Net__writeToXmlFileWithTextWriter() Description: Output a Libnucnet__Net structure to an xml file using the xml text writer. Syntax: void Libnucnet__Net__writeToXmlFileWithTextWriter( const Libnucnet__Net *self, const char *s_xml_filename, int i_indent );Input:
Upon successful return, the contents of the network are written to s_xml_filename. The routine uses the xml text writer, so it does not build a tree. Example: Dump the contents of Libnucnet__Net *p_my_network to the xml file my.xml without indentation:Libnucnet__Net__writeToXmlFileWithTextWriter( p_my_network, "my.xml", 0 ); |
Name: Libnucnet__Zone__clearNetViews() Description: Clear all the network views in a zone. Syntax: int Libnucnet__Zone__clearNetViews( Libnucnet__Zone * self );Input:
The routine returns 1 if all of the network views in the zone have been cleared. If the zone is invalid, the routine returns 0. Example: Clear the network views from p_my_zone:Libnucnet__Zone__clearNetViews( p_my_zone ); |
Name: Libnucnet__Zone__clearNseCorrectionFactorFunction() Description: Clear the NSE correction factor function for a zone. Syntax: void Libnucnet__Zone__clearNseCorrectionFactorFunction( Libnucnet__Zone *self );Input:
On return, the NSE correction factor function has been cleared and reset to the default (no correction). If the input pointer is not valid, Libnucnet error handling is invoked. Example: Clear the NSE correction factor function for Libnucnet__Zone p_zone:Libnucnet__Zone__clearNseCorrectionFactorFunction( p_zone ); |
Name: Libnucnet__Zone__clearRates() Description: Clear all the current rate values stored in a zone. Syntax: int Libnucnet__clearRates( Libnucnet__Zone * self );Input:
The routine returns 1 if all of the current rate values stored in the zone have been cleared. If the zone is invalid, the routine returns 0. Example: Clear the current rate values stored in p_my_zone:Libnucnet__clearRates( p_my_zone ); |
Name: Libnucnet__Zone__clearScreeningFunction() Description: Clear the screening function for a zone. Syntax: void Libnucnet__Zone__clearScreeningFunction( Libnucnet__Zone *self );Input:
On return, the screening function has been cleared and reset to the default (no screening). If the input pointer is not valid, Libnucnet error handling is invoked. Example: Clear the screening function for Libnucnet__Zone p_zone:Libnucnet__Zone__clearScreeningFunction( p_zone ); |
Name: Libnucnet__Zone__computeAMoment() Description: Compute the A moment (A to an integer power times abundance) for a zone. Syntax: double Libnucnet__Zone__computeAMoment( const Libnucnet__Zone *self, int n );Input:
For valid input, the routine returns the sum over all species of A of the species raised to power n times the abundance of that species in the zone. For invalid input, error handling is invoked. Example: Compute the sum of mass fractions for zone p_zone:d_xsum = Libnucnet__Zone__computeAMoment( p_zone, 1 ); |
Name: Libnucnet__Zone__computeFlowVector() Description: Creates the flow vector for the zone. Elements of the vector are the sum of all flows into the species minus the sum of all flows out of the species. The sums are computed from the individual reactions separately for best accuracy. Syntax: gsl_vector * Libnucnet__Zone__computeFlowVector( Libnucnet__Zone *self );Input:
If the input is valid, the routine returns a pointer to a new gsl_vector containing the flow vector. It is the caller's responsibility to free the vector with gsl_vector_free when no longer in use. If the input is not valid, error handling is invoked. Example: Create the flow vector for reaction rates and abundances in Libnucnet__Zone *p_zone and store it in gsl_vector *p_flow_vector:p_flow_vector = Libnucnet__Zone__computeFlowVector( p_zone ); |
Name: Libnucnet__Zone__computeJacobian() Description: Creates the Jacobian for the zone. Syntax: WnMatrix * Libnucnet__Zone__computeJacobian( Libnucnet__Zone *self );Input:
If the input is valid, the routine returns a pointer to a new WnMatrix structure containing the Jacobian. It is the caller's reponsibility to free the matrix with WnMatrix__free() when no longer in use. Example: Create the Jacobian for reaction rates and abundances in Libnucnet__Zone *p_zone and store it in WnMatrix *p_matrix:p_matrix = Libnucnet__Zone__computeJacobian( p_zone ); |
Name: Libnucnet__Zone__computeJacobianMatrix() Description: Creates the Jacobian matrix for the zone. The elements in this matrix differ from those in the regular Jacobian (computed by Libnucnet__Zone__computeJacobian) by a minus sign. Syntax: WnMatrix * Libnucnet__Zone__computeJacobianMatrix( Libnucnet__Zone *self );Input:
If the input is valid, the routine returns a pointer to a new WnMatrix structure containing the Jacobian matrix. It is the caller's reponsibility to free the matrix with WnMatrix__free() when no longer in use. Example: Create the Jacobian matrix for reaction rates and abundances in Libnucnet__Zone *p_zone and store it in WnMatrix *p_matrix:p_matrix = Libnucnet__Zone__computeJacobianMatrix( p_zone ); |
Name: Libnucnet__Zone__computeRates() Description: Compute and store the rates for all valid reactions in a given zone for the input temperature and density. Syntax: void Libnucnet__Zone__computeRates( Libnucnet__Zone *self, double d_t9, double d_rho );Input:
On successful return, the rates for each reaction stored in the input structure have been updated. Example: Update the rates for the zone with labels "0", "0", "0" for the temperature T9 = 1 and density rho = 2.e6 g/cc in Libnucnet structure p_my_zones:p_zone = Libnucnet__getZoneByLabels( p_my_zones, "0", "0", "0" ); Libnucnet__Zone__computeRates( p_zone, 1., 2.e6 ); |
Name: Libnucnet__Zone__computeZMoment() Description: Compute the Z moment (Z to an integer power times abundance) for a zone. Syntax: double Libnucnet__Zone__computeZMoment( const Libnucnet__Zone *self, int n );Input:
For valid input, the routine returns the sum over all species of Z of the species raised to power n times the abundance of that species in the zone. For invalid input, error handling is invoked. Example: Compute the Ye for zone p_zone:d_ye = Libnucnet__Zone__computeZMoment( p_zone, 1 ); |
Name: Libnucnet__Zone__copy() Description: Copy a zone. Syntax: Libnucnet__Zone * Libnucnet__Zone__copy( const Libnucnet__Zone *self );Input:
Routine returns a pointer to the new zone. If the zone could not be created, routine returns NULL. If the input zone is invalid, Libnucnet error handling is invoked. Example: Get a copy of the zone p_zone:p_copy = Libnucnet__Zone__copy( p_zone ); |
Name: Libnucnet__Zone__copy_net_views() Description: Copy the net views from the source zone to the destination zone. Syntax: void Libnucnet__Zone__copy_net_views( Libnucnet__Zone * p_destination const Libnucnet__Zone * p_source );Input:
On successful return, the network views in the source zone have been copied to the destination zone. If either zone is invalid, or if the underlying networks of the two zones are not the same, Libnucnet error handling is invoked. Example: Copy the network views from p_source to p_destination:Libnucnet__Zone__copy_net_view( p_destination, p_source ); |
Name: Libnucnet__Zone__free() Description: Free a zone. Syntax: void Libnucnet__Zone__free( Libnucnet__Zone *self );Input:
Upon successful return, the zone has been freed. If the input zone is invalid, Libnucnet error handling is invoked. Example: Free the zone p_zone:Libnucnet__Zone__free( p_zone ); |
Name: Libnucnet__Zone__getAbundanceChanges() Description: Routine to return a gsl_vector containing the changes in the abundances of a given zone. Syntax: gsl_vector * Libnucnet__Zone__getAbundanceChanges( const Libnucnet__Zone *self );Input:
Routine returns a new gsl_vector containing the current abundances in the zone. The abundances are ordered according to their index. If the input zone is not valid, Libnucnet error handling is invoked. Example: Retrieve from Libnucnet *p_my_zones the abundances in the zone with labels "x1", "y1", "z1":if( ( p_zone = Libnucnet__getZoneByLabels( p_my_zones, "x1", "y1", "z1" ) ) ) p_abund_changes = Libnucnet__Zone__getAbundanceChanges( p_zone ); |
Name: Libnucnet__Zone__getAbundances() Description: Routine to return a gsl_vector containing the abundances of a given zone. Syntax: gsl_vector * Libnucnet__Zone__getAbundances( const Libnucnet__Zone *self );Input:
Routine returns a new gsl_vector containing the current abundances in the zone. The abundances are ordered according to their index. If the input zone is not valid, Libnucnet error handling is invoked. Example: Retrieve from Libnucnet *p_my_zones the abundances in the zone with labels "x1", "y1", "z1":if( ( p_zone = Libnucnet__getZoneByLabels( p_my_zones, "x1", "y1", "z1" ) ) ) p_abunds = Libnucnet__Zone__getAbundances( p_zone ); |
Name: Libnucnet__Zone__getDataForUserRateFunction() Description: Retrieve the data for a user rate function appropriate for the given zone. Syntax: void * Libnucnet__Zone__updateDataForUserRateFunction( Libnucnet__Zone * self, const char * s_function_key );Input:
Routine returns the pointer to the user-defined data structure holding the data for the rate function. If the input zone is invalid, Libnucnet error handling is invoked. Example: Retrieve the data (originally stored in the structure my_data) for the function with key "my_function" in p_zone:p_data = ( my_data * ) Libnucnet__Zone__getDataForUserRateFunction( p_zone, "my_function" ); |
Name: Libnucnet__Zone__getEvolutionNetView() Description: Retrieve the current network view used to evolve the abundances for a zone. Syntax: Libnucnet__NetView * Libnucnet__Zone__getEvolutionNetView( Libnucnet__Zone * self );Input:
Routine returns the pointer to the Libnucnet__NetView structure with the current evolution network. If the input pointer is not valid, Libnucnet error handling is invoked. Example: Retrieve the current evolution network view for zone p_zone:p_evol_view = Libnucnet__Zone__getEvolutionNetView( p_zone ); |
Name: Libnucnet__Zone__getLabel() Description: Retrieves a label of a zone. Syntax: const char * Libnucnet__Zone__getLabel( Libnucnet__Zone *self, int i_label );Input:
Routine returns a string giving the label of the zone. If any input is invalid, Libnucnet error handling is invoked. Example: Print the second label of the zone pointed to by p_zone:printf( "Second zone label is %s\n", Libnucnet__Zone__getLabel( p_zone, 2 ) ); |
Name: Libnucnet__Zone__getMassFractions() Description: Routine to return a gsl_vector containing the mass fractions of a given zone. Syntax: gsl_vector * Libnucnet__Zone__getMassFractions( const Libnucnet__Zone *self );Input:
Routine returns a new gsl_vector containing the current mass fractions in the zone. The mass fractions are ordered according to their index. If the input zone is not valid, Libnucnet error handling is invoked. Example: Retrieve from Libnucnet *p_my_zones the mass fractions in the zone with labels "x1", "y1", "z1":if( ( p_zone = Libnucnet__getZoneByLabels( p_my_zones, "x1", "y1", "z1" ) ) ) p_x = Libnucnet__Zone__getMassFractions( p_zone ); |
Name: Libnucnet__Zone__getNet() Description: Retrieve the Libnucnet__Net structure from a Libnucnet__Zone structure. Syntax: Libnucnet__Net * Libnucnet__Zone__getNet( const Libnucnet__Zone *self );Input:
Routine returns the pointer to the Libnucnet__Net structure contained in the Libnucnet__Zone structure. If the input pointer is not valid, Libnucnet error handling is invoked. Example: Retrieve the nuclear network data structure p_my_network from the zone p_my_zone:p_my_network = Libnucnet__Zone__getNet( p_my_zone ); |
Name: Libnucnet__Zone__getNetView() Description: Retrieve the Libnucnet__NetView structure from a Libnucnet__Zone. Syntax: Libnucnet__NetView * Libnucnet__Zone__getNetView( Libnucnet__Zone * self, const char * s_label1, const char * s_label2, const char * s_label3, );Input:
Routine returns the pointer to the Libnucnet__NetView structure. If the input pointer is not valid, Libnucnet error handling is invoked. Example: Retrieve the network view stored with label "my view" from Zone p_zone:p_my_view = Libnucnet__getNetView( p_zone, "my view", NULL, NULL ); |
Name: Libnucnet__Zone__getNseCorrectionFactorData() Description: Retrieve the user-supplied NSE correction factor function data for a zone. Syntax: void * Libnucnet__Zone__getNseCorrectionFactorData( Libnucnet__Zone *self );Input:
For valid input, the routine returns the data for the NSE correction factor function for the zone. Example: Get the user-supplied data for the zone p_zone:p_data = Libnucnet__Zone__getNseCorrectionFactorData( p_zone ); |
Name: Libnucnet__Zone__getNseCorrectionFactorFunction() Description: Retrieve the user-supplied NSE correction factor function for a zone. Syntax: Libnucnet__Species__nseCorrectionFactorFunction Libnucnet__Zone__getNseCorrectionFactorFunction( Libnucnet__Zone *self );Input:
For valid input, the routine returns the NSE correction factor function for the zone. Example: Get the user-supplied NSE correction factor routine for the zone p_zone:Libnucnet__Species__nseCorrectionFactorFunction pf_correction_function = Libnucnet__Zone__getNseCorrectionFactorFunction( p_zone ); |
Name: Libnucnet__Zone__getProperty() Description: Retrieve a zone optional property. Syntax: const char * Libnucnet__Zone__getProperty( const Libnucnet__Zone *self, const char *s_name, const char *s_tag1, const char *s_tag2 );Input:
Routine returns a const char * containing the value of the property. If the property is not found, routine returns NULL. Example: Print the optional property t9 in p_zone:printf( "t9 = %g\n", atof( Libnucnet__Zone__getProperty( p_zone, "t9", NULL, NULL ) ) ); |
Name: Libnucnet__Zone__getRatesForReaction() Description: Retrieve the rates for the specified reaction from the input zone. Syntax: void Libnucnet__Zone__getRatesForReaction( Libnucnet__Zone *self, Libnucnet__Reaction *p_reaction, double *p_forward, double *p_reverse );Input:
p_zone = Libnucnet__getZoneByLabels( p_zones, "1", "1", "2" ); Libnucnet__Zone__getRatesForReaction( p_zone, p_reaction, p_forward, p_reverse ); |
Name: Libnucnet__Zone__getScreeningData() Description: Retrieve the user-supplied screening factor function data for a zone. Syntax: void * Libnucnet__Zone__getScreeningData( Libnucnet__Zone *self );Input:
For valid input, the routine returns the data for the screening function for the zone. Example: Get the user-supplied screening routine for the zone p_zone:p_data = Libnucnet__Zone__getScreeningData( p_zone ); |
Name: Libnucnet__Zone__getScreeningFunction() Description: Retrieve the user-supplied screening factor function for a zone. Syntax: Libnucnet__Zone__screeningFunction Libnucnet__Zone__getScreeningFunction( Libnucnet__Zone *self );Input:
For valid input, the routine returns the screening function for the zone. Example: Get the user-supplied screening routine for the zone p_zone:Libnucnet__Zone__screeningFunction pf_screening_function = Libnucnet__Zone__getScreeningFunction( p_zone ); |
Name: Libnucnet__Zone__getSpeciesAbundance() Description: Retrieves the specified species abundance. Syntax: double Libnucnet__Zone__getSpeciesAbundance( const Libnucnet__Zone *self, const Libnucnet__Species *p_species );Input:
Routine returns a double giving the abundance of the species in the specified zone. If any input is invalid, Libnucnet error handling is invoked. Example: Retrieve he4 abundance in zone p_zone from the network p_my_nucnet and store in d_abundance:d_abundance = Libnucnet__Zone__getSpeciesAbundance( p_zone, Libnucnet__Nuc__getSpeciesByName( Libnucnet__Net__getNuc( Libnucnet__Zone__getNet( p_my_nucnet ) ), "he4" ) ); |
Name: Libnucnet__Zone__getSpeciesAbundanceChange() Description: Retrieves the specified species abundance change. Syntax: double Libnucnet__Zone__getSpeciesAbundanceChange( const Libnucnet__Zone *self, const Libnucnet__Species *p_species );Input:
Routine returns a double giving the change in abundance of the species in the specified zone. If any input is invalid, Libnucnet error handling is invoked. Example: Retrieve he4 abundance change in zone p_zone:if( ( p_he4 = Libnucnet__Nuc__getSpeciesByName( Libnucnet__Net__getNuc( Libnucnet__Zone__getNet( p_zone ) ), "he4" ) ) ) { d_abundance_change = Libnucnet__Zone__getSpeciesAbundanceChange( p_zone, p_he4 ); } |
Name: Libnucnet__Zone__getSpeciesMassFraction() Description: Retrieves the specified species mass fraction. Syntax: double Libnucnet__Zone__getSpeciesMassFraction( const Libnucnet__Zone *self, const Libnucnet__Species *p_species );Input:
Routine returns a double giving the mass fraction of the species in the specified zone. If any input is invalid, Libnucnet error handling is invoked. Example: Retrieve he4 mass fraction in zone p_zone from the network p_my_nucnet and store in d_x:d_x = Libnucnet__Zone__getSpeciesMassFraction( p_zone, Libnucnet__Nuc__getSpeciesByName( Libnucnet__Net__getNuc( Libnucnet__Zone__getNet( p_my_nucnet ) ), "he4" ) ); |
Name: Libnucnet__Zone__getSummedAbundances() Description: Retrieve the abundances for a zone as a function of Z, N, or A. Syntax: gsl_vector * Libnucnet__Zone__getSummedAbundances( Libnucnet__Zone * self, const char * s_nucleon_type );Input:
Routine returns a new gsl_vector that gives the abundances in the zone summed by Z, N, or A. The index of a component of the vector is the particular value of Z, N, or A while the value of the component is the summed abundance over that index. If the input zone or input string for the nucleon number type is invalid, error handling is invoked. Example: Get the abundances for Libnucnet__Zone * p_zone as a function of atomic number:p_z_vector = Libnucnet__Zone__getSummedAbundances( p_zone, "z" ); |
Name: Libnucnet__Zone__isComputingReverseRatesFromDetailedBalance() Description: Determine whether a zone is computing reverse reaction rates from detailed balance on the forward rates. Syntax: int Libnucnet__Zone__isComputingReverseRatesFromDetailedBalance( const Libnucnet__Zone * self );Input:
Routine returns 1 (true) if the zone is computing reverse rates from detailed balance and 0 (false) if not. If the input zone is invalid, error handling is invoked. Example: Determine whether the zone p_zone is computing reverse rates from detailed balance:if( Libnucnet__Zone__isComputingReverseRatesFromDetailedBalance( p_zone ) ) fprintf( stdout, "Zone is using detailed balance for reverse rates.\n" ); |
Name: Libnucnet__Zone__iterateNetViews() Description: Iterate over the network views in a zone. Syntax: void Libnucnet__Zone__iterateNetViews( const Libnucnet__Zone *self, const char *s_label1, const char *s_label2, const char *s_label3, Libnucnet__Zone__NetView_iterateFunction pf_func, void *p_user_data );Input:
Routine iterates over the network views in the zone and applies the user-supplied function. The properties iterated over are those that match s_label1, s_label2, and s_label3. If any of these is NULL, the comparison is a match. The views are iterated in the order they are stored internally over which the user has no control. Examples: Iterate over all the network views in zone p_zone and apply my_func and the extra data pointed to by p_user_data:Libnucnet__Zone__iterateNetViews( p_zone, NULL, NULL, NULL, (Libnucnet__NetView__iterateFunction) my_func, p_user_data );Iterate over all the network views that match "my views" as the first label in zone p_zone and apply my_func and the extra data pointed to by p_user_data: Libnucnet__Zone__iterateNetViews( p_zone, "my views", NULL, NULL, (Libnucnet__NetView__iterateFunction) my_func, p_user_data ); |
Name: Libnucnet__Zone__iterateOptionalProperties() Description: Iterate over the optional properties in a zone. Syntax: void Libnucnet__Zone__iterateOptionalProperties( const Libnucnet__Zone *self, const char *s_name, const char *s_tag1, const char *s_tag2, Libnucnet__Zone__optional_property_iterate_function pf_func, void *p_user_data );Input:
Routine iterates over the optional properties in the zone and applies the user-supplied function. The properties iterated over are those that match s_name, s_tag1, and s_tag2. If any of these is NULL, the comparison is a match. The properties are iterated in the order they are stored internally over which the user has no control. Examples: Iterate over the all optional properties in zone p_zone and apply my_func and the extra data pointed to by p_user_data:Libnucnet__Zone__iterateOptionalProperties( p_zone, NULL, NULL, NULL, (Libnucnet__Zone__optional_property_iterate_function) my_func, p_user_data );Iterate over all the optional properties that match "t9" as a name in zone p_zone and apply my_func and the extra data pointed to by p_user_data: Libnucnet__Zone__iterateOptionalProperties( p_zone, "t9", NULL, NULL, (Libnucnet__Zone__optional_property_iterate_function) my_func, p_user_data ); |
Name: Libnucnet__Zone__new() Description: Create a new zone. Syntax: Libnucnet__Zone * Libnucnet__Zone__new( Libnucnet__Net *p_net, const char *s_label1, const char *s_label2, const char *s_label3 );Input:
Routine returns a pointer to the new zone. If the zone could not be created, routine returns NULL. If any other input is invalid, Libnucnet error handling is invoked. Example: Create a new 1 dimensional zone with label "Z" with data from the network p_net:p_new_zone = Libnucnet__Zone__new( p_net, "Z", NULL, NULL ); |
Name: Libnucnet__Zone__removeNetView() Description: Remove a network view from a zone. Syntax: int Libnucnet__Zone__removeNetView( Libnucnet__Zone * self, const char * s_label1, const char * s_label2, const char * s_label3 );Input:
Returns 1 (true) if the removal was successful and 0 (false) if not. On successful return, the network view has been removed from the zone and the memory for the view has been cleared. If input is invalid, Libnucnet error handling is invoked. Example: Remove the network view labeled "my view" from p_zone:if( Libnucnet__Zone__removeNetView( p_zone, "my view", NULL, NULL ) ) fprintf( stdout, "Removal succeeded.\n" ); |
Name: Libnucnet__Zone__removeProperty() Description: Remove a zone optional property. Syntax: int Libnucnet__Zone__removeProperty( Libnucnet__Zone *self, const char *s_name, const char *s_tag1, const char *s_tag2 );Input:
Routine returns 1 if removal was successful and 0 if not. On successful return, the property has been removed from the zone. If the input are invalid, error handling is invoked. Example: Remove the optional property t9 from p_zone:if !Libnucnet__Zone__removeProperty( p_zone, "t9", NULL, NULL, ) ) fprintf( stderr, "Unsuccessful removal!\n" ); |
Name: Libnucnet__Zone__setNseCorrectionFactorFunction() Description: Set the user-supplied correction factor function and its associated data. Syntax: void Libnucnet__Zone__setNseCorrectionFactorFunction( Libnucnet__Zone *self, Libnucnet__Species__nseCorrectionFactorFunction pf_user_function, void *p_user_data );Input:
For valid input, the routine sets the correction factor function to that supplied and the data to the user's data structure. These are applied in the routine Libnucnet__Zone__computeRates(). Examples: Set the user-supplied correction factor routine to my_correction_function and pass no extra data for rate calculations in Libnucnet__Zone p_zone:Libnucnet__Zone__setNseCorrectionFactorFunction( p_zone, (Libnucnet__Species__nseCorrectionFactorFunction) my_correction_function, NULL );Set the user-supplied correction factor routine to my_correction_function and pass extra data from the structure my_data for rate calculations in Libnucnet__Zone p_zone: struct my_data = { double x; double y; }; my_data.x = 1.; my_data.y = 2.; Libnucnet__Zone__setNseCorrectionFactorFunction( p_zone, (Libnucnet__Species__nseCorrectionFactorFunction) my_correction_function, &my_data ); |
Name: Libnucnet__Zone__setScreeningFunction() Description: Set the user-supplied screening factor function and its associated data. Syntax: void Libnucnet__Zone__setScreeningFunction( Libnucnet__Zone *self, Libnucnet__Zone__screeningFunction pf_user_function, void *p_user_data );Input:
For valid input, the routine sets the screening function to that supplied and the data to the user's data structure. These are applied in the routine Libnucnet__Zone__computeRates(). Examples: Set the user-supplied screening routine to my_screening_function and pass no extra data for rate calculations in Libnucnet__Zone p_zone:Libnucnet__Zone__setScreeningFunction( p_zone, (Libnucnet__Zone__screeningFunction) my_screening_function, NULL );Set the user-supplied screening routine to my_screening_function and pass extra data from the structure my_data for rate calculations in Libnucnet__Zone p_zone: struct my_data = { double x; double y; }; my_data.x = 1.; my_data.y = 2.; Libnucnet__Zone__setScreeningFunction( p_zone, (Libnucnet__Zone__screeningFunction) my_screening_function, &my_data ); |
Name: Libnucnet__Zone__toggleReverseRateDetailedBalance() Description: Turn on or off computation of the reverse rates for reactions from detailed balance for a zone. Syntax: void Libnucnet__Zone__toggleReverseRateDetailedBalance( Libnucnet__Zone * self, const char * s_switch );Input:
On successful return, the zone will be set to compute reverse rates by detailed balance or not. If the input zone is not valid, or if s_switch is not one of "on" or "off", Libnucnet error handling is invoked. The default setting for a zone is "on". Example: Turn off computation of the reverse rates for p_zone:Libnucnet__Zone__toggleReverseRateDetailedBalance( p_zone, "off" ); |
Name: Libnucnet__Zone__updateAbundanceChanges() Description: Routine to update the abundance changes of a given zone from an input gsl_vector. Syntax: void Libnucnet__Zone__updateAbundanceChanges( Libnucnet__Zone *self gsl_vector *p_abunds );Input:
Upon successful return, the abundance changes have been updated to those in the input vector. If the input zone is not valid, or if the number of elements in the input vector is not equal to the number of species in the zone's network, Libnucnet error handling is invoked. Example: Update the abundance changes in Libnucnet__Zone *p_zone with the abundance changes in gsl_vector *p_my_new_abundance_changes:Libnucnet__Zone__updateAbundanceChanges( p_zone, p_my_new_abundance_changes ); |
Name: Libnucnet__Zone__updateAbundances() Description: Routine to update the abundances of a given zone from an input gsl_vector. Syntax: void Libnucnet__Zone__updateAbundances( Libnucnet__Zone *self gsl_vector *p_abunds );Input:
Upon successful return, the abundances have been updated to those in the input vector. If the input zone is not valid, or if the number of elements in the input vector is not equal to the number of species in the zone's network, Libnucnet error handling is invoked. Example: Update the abundances in Libnucnet__Zone *p_zone with the abundances in gsl_vector *p_my_new_abundances:Libnucnet__Zone__updateAbundances( p_zone, p_my_new_abundances ); |
Name: Libnucnet__Zone__updateDataForUserRateFunction() Description: Update the data for a user rate function appropriate for the given zone. Syntax: int Libnucnet__Zone__updateDataForUserRateFunction( Libnucnet__Zone * self, const char * s_function_key, void * p_data );Input:
Routine returns 1 (true) if the update succeeded and 0 (false) if not. If the input zone is invalid, Libnucnet error handling is invoked. Example: Update the function with key "my_function" with data in structure my_data for zone p_zone:Libnucnet__Zone__updateDataForUserRateFunction( p_zone, "my_function", &my_data ); |
Name: Libnucnet__Zone__updateMassFractions() Description: Routine to update the mass fractions of a given zone from an input gsl_vector. Syntax: void Libnucnet__Zone__updateMassFractions( Libnucnet__Zone *self gsl_vector *p_x );Input:
Upon successful return, the mass fractions have been updated to those in the input vector. If the input zone is not valid, or if the number of elements in the input vector is not equal to the number of species in the zone's network, Libnucnet error handling is invoked. Example: Update the abundances in Libnucnet__Zone *p_zone with the mass fractions in gsl_vector *p_my_new_mass_fractions:Libnucnet__Zone__updateMassFractions( p_zone, p_my_new_mass_fractions ); |
Name: Libnucnet__Zone__updateNetView() Description: Update a network view in a zone. Syntax: int Libnucnet__Zone__updateNetView( Libnucnet__Zone * self, const char * s_label1, const char * s_label2, const char * s_label3, Libnucnet__NetView * p_view );Input:
Returns 1 (true) if the update was successful and 0 (false) if not. On successful return, the network view has been added to the zone if the zone previously did not include the view or has been updated in the zone. If input is invalid, Libnucnet error handling is invoked. Examples: Update the network view labeled "my nuclides" and "my reactions" in p_zone with view p_view:if( Libnucnet__Zone__updateNetView( p_zone, "my nuclides", "my reactions", NULL, p_view ) ) fprintf( stdout, "Update succeeded.\n" );Update the evolution network view for p_zone with view p_view: if( Libnucnet__Zone__updateNetView( p_zone, EVOLUTION_NETWORK, NULL, NULL, p_view ) ) fprintf( stdout, "Update succeeded.\n" ); |
Name: Libnucnet__Zone__updateProperty() Description: Update a zone optional property. Syntax: int Libnucnet__Zone__updateProperty( Libnucnet__Zone *self, const char *s_name, const char *s_tag1, const char *s_tag2, const char *s_value );Input:
Routine returns 1 if the update was successful and 0 if not. On successful return, the property has been added to the zone if it did not previously exist or has been updated if it did. If the input are invalid, error handling is invoked. Example: Update the optional property t9 to 4.5 in p_zone:Libnucnet__Zone__updateProperty( p_zone, "t9", NULL, NULL, "4.5" ); |
Name: Libnucnet__Zone__updateRatesForReaction() Description: Updates the rates for a given reaction to the specified values in the input zone. Syntax: void Libnucnet__Zone__updateRatesForReaction( Libnucnet__Zone *self, Libnucnet__Reaction *p_reaction, double d_forward, double d_reverse );Input:
p_zone = Libnucnet__getZoneByLabels( p_zones, "1", "2", "3" ); Libnucnet__Zone__getRatesForReaction( p_zone, p_reaction, &d_forward, &d_reverse ); Libnucnet__Zone__updateRatesForReaction( p_zone, p_reaction, 2. * d_forward, 2. * d_reverse ); |
Name: Libnucnet__Zone__updateSpeciesAbundance() Description: Assigns the new species abundance to the Libnucnet__Zone struct. Syntax: void Libnucnet__Zone__updateSpeciesAbundance( Libnucnet__Zone *self, Libnucnet__Species * p_species, double d_y );Input:
On successful return, the abundance of the species in the zone is modified by the value in d_y. Example: Update the abundance of species p_species for p_zone with d_y:Libnucnet__Zone__updateSpeciesAbundance( p_zone, p_species, d_y ); |
Name: Libnucnet__Zone__updateSpeciesAbundanceChange() Description: Assigns the new species abundance change to the Libnucnet__Zone struct. Syntax: void Libnucnet__Zone__updateSpeciesChangeAbundance( Libnucnet__Zone *self, Libnucnet__Species * p_species, double d_dy );Input:
On successful return, the abundance change of the species in the zone is modified by the value in d_dy. Example: Update the abundance change of species p_species for p_zone with d_dy:Libnucnet__Zone__updateSpeciesAbundanceChange( p_zone, p_species, d_dy ); |
Name: Libnucnet__Zone__updateSpeciesMassFraction() Description: Assigns the new species mass fraction to the Libnucnet__Zone struct. Syntax: void Libnucnet__Zone__updateSpeciesMassFraction( Libnucnet__Zone *self, Libnucnet__Species * p_species, double d_x );Input:
On successful return, the mass fraction of the species in the zone is modified by the value in d_x. Example: Update the mass fraction of species p_species for p_zone with d_x:Libnucnet__Zone__updateSpeciesMassFraction( p_zone, p_species, d_x ); |
Name: Libnucnet__Zone__updateTimeStep() Description: Calculates the new time step for the evolution of the nuclear network. Syntax: void Libnucnet__Zone__updateTimeStep( Libnucnet__Zone *self, double *d_dt, double d_regt, double d_regy, double d_ymin );Input:
On successful return, the timestep has been updated. If input is not valid, error handling is invoked. Example: Update the time step for p_zone given old timestep of d_dt, a maximum allowed increase in the timestep of 15%, and a maximum allowed increase in the abundance of a species with abundance greater than 1.e-10 of 10%:Libnucnet__Zone__updateTimeStep( p_zone, &d_dt, 0.15, 0.1, 1.e-10 ); |
Name: Libnucnet__addZone() Description: Add a zone to a Libnucnet structure. Syntax: int Libnucnet__addZone( Libnucnet *self, Libnucnet__Zone *p_zone );Input:
Routine returns 1 if the addition succeeded and 0 if not. If any other input is invalid, Libnucnet error handling is invoked. Example: Add p_zone to Libnucnet *p_my_nucnet:if( !Libnucnet__addZone( p_my_nucnet, p_zone ) ) printf( "Zone addition failed!\n" ); |
Name: Libnucnet__assignZoneDataFromXml() Description: Reads in zone data from an XML file, creates zones, and stores the data in a Libnucnet structure. Syntax: void Libnucnet__assignZoneDataFromXml( Libnucnet *self, const char *s_xml_filename, const char *s_zone_xpath );Input:
If the input is valid, the Libnucnet structure has been updated with the zones. Examples: Store the initial mass fraction data in input_mass.xml in p_my_nucnet:Libnucnet__assignZoneDataFromXml( p_my_nucnet, "input_mass.xml", NULL );Store the initial mass fraction data from mass_fractions.xml but only from zones with value 0 for label 1: Libnucnet__assignZoneDataFromXml( p_my_nucnet, "mass_fractions.xml", "[@label1='0']" ); |
Name: Libnucnet__assignZoneDataFromXmlTextReader() Description: Reads in zone data from an XML file, creates zones, and stores the data in a Libnucnet structure without building a tree. Syntax: void Libnucnet__assignZoneDataFromXmlTextReader( Libnucnet *self, const char *s_xml_filename, const char *s_zone_xpath );Input:
If the input is valid, the Libnucnet structure has been updated with the zones. Examples: Store the initial mass fraction data in input_mass.xml in p_my_nucnet:Libnucnet__assignZoneDataFromXmlTextReader( p_my_nucnet, "input_mass.xml", NULL );Store the initial mass fraction data from mass_fractions.xml but only from zones with value 0 for label 1: Libnucnet__assignZoneDataFromXml( p_my_nucnet, "mass_fractions.xml", ".//@label1[.='0']" ); |
Name: Libnucnet__clearZoneCompareFunction() Description: Restore the default function to be applied during a zone sort. Syntax: void Libnucnet__clearZoneCompareFunction( Libnucnet *self );Input:
Upon successful return, the zone compare function for the collection has been restored to the default function. If any input is invalid, error handling is invoked. Example: Clear the zone compare function in p_my_nucnet:Libnucnet__clearZoneCompareFunction( p_my_nucnet ); |
Name: Libnucnet__free() Description: Free the memory for a Libnucnet structure. Syntax: void Libnucnet__free( Libnucnet *self );Input:
On return, the memory allocated for the network has been cleared. If the input pointer is not valid, Libnucnet error handling is invoked. Example: Free the memory for the network in p_my_nucnet:Libnucnet__free( p_my_nucnet ); |
Name: Libnucnet__freeAllZones() Description: Free the memory for all zones in a Libnucnet structure. Syntax: void Libnucnet__freeAllZones( Libnucnet *self );Input:
On return, the memory allocated for all zones has been cleared but self is still available for use. If the input pointer is not valid, Libnucnet error handling is invoked. Example: Free the memory for all zones in p_my_nucnet:Libnucnet__freeAllZones( p_my_nucnet ); |
Name: Libnucnet__getNet() Description: Retrieve the Libnucnet__Net structure from a Libnucnet structure. Syntax: Libnucnet__Net *Libnucnet__getNet( Libnucnet *self );Input:
Routine returns the pointer to the Libnucnet__Net structure contained in the Libnucnet structure. If the input pointer is not valid, Libnucnet error handling is invoked. Example: Retrieve the nuclear network data structure p_my_network from the libnucnet structure p_my_nucnet:p_my_network = Libnucnet__getNet( p_my_nucnet ); |
Name: Libnucnet__getNumberOfZones() Description: Retrieves the number of zones in the network. Syntax: size_t Libnucnet__getNumberOfZones( const Libnucnet *self );Input:
Routine returns an integer representing the number of zones in the network. If the input stucture is not valid, error handling is invoked. Example: Print number of zones in Libnucnet *p_my_nucnet:printf( "Number of zones = %d\n", Libnucnet__getNumberOfZones( p_my_nucnet ) ); |
Name: Libnucnet__getZoneByLabels() Description: Retrieves the specified zone. Syntax: Libnucnet__Zone * Libnucnet__getZoneByLabels( const Libnucnet *self, const char *s_label_1, const char *s_label_2, const char *s_label_3 );Input:
Routine returns a pointer to a Libnucnet__Zone structure. If the input is invalid, Libnucnet error handling is invoked. If the zone is not found, routine returns NULL. Example: Retrieve zone with labels "x", "y", and "z" from Libnucnet *p_nucnet:p_zone = Libnucnet__getZoneByLabels( p_nucnet, "x", "y", "z" ); |
Name: Libnucnet__is_valid_input_xml() Description: Validate an xml file for Libnucnet. Syntax: int Libnucnet__is_valid_input_xml( const char *s_xml_filename );Input:
For a valid input Libnucnet data xml file, the routine returns 1 (true). For an invalid file, routine returns 0 (false). If the schema file is invalid, or if the schema file cannot be read over the web, routine stops and prints error message. Example: Validate the input xml file "libnucnet_input.xml":if( Libnucnet__is_valid_input_xml( "libnucnet_input.xml" ) ) { printf( "Valid xml!\n" ); } |
Name: Libnucnet__is_valid_zone_data_xml() Description: Validate an xml file for zone data. Syntax: int Libnucnet__is_valid_zone_data_xml( const char *s_xml_filename );Input:
For a valid input zone data xml file, the routine returns 1 (true). For an invalid file, routine returns 0 (false). If the schema file is invalid, or if the schema file cannot be read over the web, routine stops and prints error message. Example: Validate the input xml file "my_zone_data.xml":if( Libnucnet__is_valid_zone_data_xml( "my_zone_data.xml" ) ) { printf( "Valid xml!\n" ); } |
Name: Libnucnet__iterateZones() Description: Iterate through the zones and apply the user-supplied iterate function. Syntax: void Libnucnet__iterateZones( const Libnucnet *self, Libnucnet__Zone__iterateFunction pf_func, void *p_user_data );Input:
The routine iterates through the zones and applies the user-supplied routine to each. If any input is invalid, error handling is invoked. Example: Iterate through the zones in p_my_nucnet and apply the function my_iterate_function and the extra data in p_user_data:Libnucnet__iterateZones( p_my_nucnet, (Libnucnet__Zone__iterateFunction) my_iterate_function, p_user_data ); |
Name: Libnucnet__new() Description: Creates a new Libnucnet structure. Syntax: Libnucnet *Libnucnet__new( );Output: Routine returns a pointer to a Libnucnet structure. If the routine cannot allocate sufficient memory, Libnucnet error handling is invoked. Example: Create a new network structure p_my_nucnet:p_my_nucnet = Libnucnet__new( ); |
Name: Libnucnet__new_from_xml() Description: Reads in nuclear and reaction data from xml file and stores the reactions in a Libnucnet structure. The reactions to be stored may be selected with an xpath expression. The routine also checks that all reactants and products in a reaction are included in the network and flags that reaction if not. Syntax: Libnucnet * Libnucnet__new_from_xml( const char *s_xml_filename, char *s_nuc_xpath char *s_reac_xpath char *s_zone_xpath );Input:
If the input is valid, routine returns a pointer to a Libnucnet structure containing nuclides, reactions, and zones selected by the xpath expressions. If any of the xpath expressions is invalid, error handling is invoked. Examples: Store all the data from "input_data.xml" to p_my_nucnet:p_my_nucnet = Libnucnet__new_from_xml( "input_data.xml", NULL, NULL, NULL );Store the data from "input_data.xml", but only those reactions, containing 'h1' as a reactant, to p_my_nucnet: p_my_nucnet = Libnucnet__new_from_xml( "input_data.xml", NULL, "[reactant = 'h1']", NULL ); |
Name: Libnucnet__new_from_xml_text_reader() Description: Reads in nuclear and reaction data from xml file and stores the reactions in a Libnucnet structure without building a tree. The reactions to be stored may be selected with an xpath expression. The routine also checks that all reactants and products in a reaction are included in the network and flags that reaction if not. For this routine, xpath expressions are relative to individual nuclides, reactions, or zones. Syntax: Libnucnet * Libnucnet__new_from_xml_text_reader( const char *s_xml_filename, char *s_nuc_xpath char *s_reac_xpath char *s_zone_xpath );Input:
If the input is valid, routine returns a pointer to a Libnucnet structure containing nuclides, reactions, and zones selected by the xpath expressions. If any of the xpath expressions is invalid, error handling is invoked. Examples: Store all the data from "input_data.xml" to p_my_nucnet:p_my_nucnet = Libnucnet__new_from_xml( "input_data.xml", NULL, NULL, NULL );Store the data from "input_data.xml", but only those reactions, containing 'h1' as a reactant, to p_my_nucnet: p_my_nucnet = Libnucnet__new_from_xml( "input_data.xml", NULL, ".//reactant[. = 'h1']", NULL ); |
Name: Libnucnet__relabelZone() Description: Change the labels of a zone. Syntax: int Libnucnet__relabelZone( Libnucnet *self, Libnucnet__Zone *p_zone, const char *s_new_label1, const char *s_new_label2, const char *s_new_label3 );Input:
Routine returns 1 if the relabeling was successful and 0 if not. On successful return, the zone has been relabeled. On unsuccessful return (for example, if a zone with the new labels already exists in self), the labels have not changed. If self or p_zone is invalid, error handling is invoked. Example: Change the labels of a zone in p_my_nucnet from "x1", "y1", and "z1" to "x2", "y2", and "z2":p_zone = Libnucnet__getZoneByLabels( p_my_nucnet, "x1", "y1", "z1" ); if( Libnucnet__relabelZone( p_my_nucnet, p_zone, "x2", "y2"," z2" ) ) printf( "Relabeling sucessful!\n" ); |
Name: Libnucnet__removeZone() Description: Remove a zone from a Libnucnet structure. Syntax: void Libnucnet__removeZone( Libnucnet *self, Libnucnet__Zone *p_zone );Input:
Upon successful return, the zone has been removed from the Libnucnet structure. If the input zone is invalid, Libnucnet error handling is invoked. Example: Remove the zone labeled "x1", "y1", "z1" from Libnucnet *p_my_nucnet:Libnucnet__removeZone( p_my_nucnet, Libnucnet__getZoneByLabels( p_my_nucnet, "x1", "y1", "z1" ) ); |
Name: Libnucnet__setZoneCompareFunction() Description: Set the function to be applied to sort the zones during a zone iteration. Syntax: void Libnucnet__setZoneCompareFunction( Libnucnet *self, Libnucnet__Zone__compare_function pfFunc );Input:
Upon successful return, the data compare function for zones in the collection has been set to the input function. If any input is invalid, error handling is invoked. Example: Set the zone compare function in p_my_nucnet to my_sort_function:Libnucnet__setZoneCompareFunction( p_my_nucnet, (Libnucnet__Zone__compare_function) my_sort_function ); |
Name: Libnucnet__updateZoneXmlMassFractionFormat() Description: Update the format code for output of zone mass fractions to xml. Syntax: void Libnucnet__updateZoneXmlMassFractionFormat( Libnucnet * self, const char * s_format );Input:
On successful return, the format code for the xml output of the zone mass fractions has been updated. If the input structure is invalid, Libnucnet error handling is invoked. Example: Update the format code for Libnucnet *p_my_nucnet for xml output of zone mass fractions to %.15e:Libnucnet__updateZoneXmlMassFractionFormat( p_my_nucnet, "%.15e" ); |
Name: Libnucnet__writeToXmlFile() Description: Output a Libnucnet structure to an xml file. Syntax: void Libnucnet__writeToXmlFile( const Libnucnet *self, const char *s_xml_filename );Input:
Upon successful return, the contents of the full libnucnet structure are written to s_xml_filename. Example: Dump the contents of Libnucnet *p_my_nucnet to the xml file my.xml:Libnucnet__writeToXmlFile( p_my_nucnet, "my.xml" ); |
Name: Libnucnet__writeToXmlFileWithTextWriter() Description: Output a Libnucnet structure to an xml file. The routine does not build an internal xml tree. Syntax: void Libnucnet__writeToXmlFileWithTextWriter( const Libnucnet *self, const char *s_xml_filename, int indent );Input:
Upon successful return, the contents of the full libnucnet structure are written to s_xml_filename. Example: Dump the contents of Libnucnet *p_my_nucnet to the xml file my.xml with indentation:Libnucnet__writeToXmlFileWithTextWriter( p_my_nucnet, "my.xml", 1 ); |
Name: Libnucnet__writeZoneDataToXmlFile() Description: Dump the current data in all zones to an xml file. Syntax: void Libnucnet__writeZoneDataToXmlFile( const Libnucnet *self const char *s_xml_file );Input:
Upon successful return, the mass fractions have been dumped to the xml file. Zones will be written out in the order defined by the currently defined Libnucnet__Zone__compare_function. If the input structure is not valid, or the xml file cannot be created, error handling is invoked. Example: Dump the current zone data in p_my_libnucnet to the file zone_data.xml.Libnucnet__writeZoneDataToXmlFile( p_my_libnucnet, "zone_data.xml" ); |
Name: Libnucnet__writeZoneDataToXmlFileWithTextWriter() Description: Dump the current data in all zones to an xml file. Don't build a tree. Syntax: void Libnucnet__writeZoneDataToXmlFileWithTextWriter( const Libnucnet *self const char *s_xml_file, int i_indent );Input:
Upon successful return, the mass fractions have been dumped to the xml file. Zones will be written out in the order defined by the currently defined Libnucnet__Zone__compare_function. If the input structure is not valid, or the xml file cannot be created, error handling is invoked. The routine does not build a true, so memory overhead is low. Example: Dump the current zone data in p_my_libnucnet to the file zone_data.xml with indentation.Libnucnet__writeZoneDataToXmlFileWithTextWriter( p_my_libnucnet, "zone_data.xml", 1 ); |