Warning: This is no longer the latest available version of this module. Please see the releases page for the most recent version. The Webnucleo group strongly recommends the use of the latest version of any of its online modules.
Documentation from the file
Table of Contents
-
Structures
-
User-Supplied Routines
-
Routines
| Structures |
|---|
Name: Libnucnet__NetDescription: Libnucnet__Net is a structure that stores data about the nuclei and reactions among them in a nuclear reaction network. A Libnucnet__Nuc structure and a Libnucnet__Reac structure comprise the Libnucnet__Net structure. The contents of Libnucnet__Net are not made public by the API but rather are accessed through the API routines. |
|
|
Name: Libnucnet__ZoneDescription: 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: LibnucnetDescription: 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__Net__screeningFunction()Description: Optional user-supplied screening function. Syntax:
double
Libnucnet__Net__screeningFunction(
double d_t9,
double d_rho,
double d_ye,
int i_Z1,
int i_A1,
int i_Z2,
int i_A2,
void *p_user_data
);
Input:
User's routine must return a double giving the screening factor for the reaction (i_Z1,i_A1) + (i_Z2,i_A2) at the input temperature, density, and Ye. The user can call this routine by passing it into Libnucnet__Net__computeScreeningFactorForReaction() or can set it with Libnucnet__Zone__setScreeningFunction for libnucnet to use from Libnucnet__Zone__computeRates(). |
|
|
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. The user's routine must not modify self, that is, the pointer to the input Libnucnet__Zone structure. |
|
|
| Routines |
|---|
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,
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,
&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__computeScreeningFactorForReaction()Description: Use the user-supplied screening factor function and its associated data to compute the screening factor for a reaction. Syntax:
double
Libnucnet__Net__computeScreeningFactorForReaction(
Libnucnet__Net *self,
Libnucnet__Reaction *p_reaction,
double d_t9,
double d_rho,
double d_ye,
Libnucnet__Net__screeningFunction p_user_function,
void *p_user_data
);
Input:
For valid input, the routine returns a double with the screening factor as computed by the user's supplied function and data. Examples: Compute the screening 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_screening_function and pass no extra data for the screening calculation:
p_reaction =
Libnucnet__Reac__getReactionByString(
Libnucnet__Net__getReac( p_my_net ),
"c12 + h1 -> n13 + gamma"
);
d_screen =
Libnucnet__Net__computeScreeningFactorForReaction(
p_my_net,
p_reaction,
2.,
1.e11,
0.5,
(Libnucnet__Net__screeningFunction) my_screening_function,
NULL
);
Compute the screening 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_screening_function
and pass extra data for the screening calculation through the
data structure my_data:
struct my_data = { double x = 1.; double y = 2.; };
p_reaction =
Libnucnet__Reac__getReactionByString(
Libnucnet__Net__getReac( p_my_net ),
"c12 + h1 -> n13 + gamma"
);
d_screen =
Libnucnet__Net__computeScreeningFactorForReaction(
p_my_net,
p_reaction,
2.,
1.e11,
0.5,
(Libnucnet__Net__screeningFunction) my_screening_function,
&my_data
);
|
|
|
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( 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( 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(
Libnucnet__Net *self,
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__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__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(
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 reactions 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__computeJacobianMatrix()Description: Creates the Jacobian matrix for the zone. 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(
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__getAbundanceChanges()Description: Routine to return a gsl_vector containing the changes in the abundances of a given zone. Syntax:
gsl_vector *
Libnucnet__Zone__getAbundanceChanges( 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( 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__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__getNet()Description: Retrieve the Libnucnet__Net structure from a Libnucnet__Zone structure. Syntax:
Libnucnet__Net *
Libnucnet__Zone__getNet( 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__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__getReverseRatioCorrectionFactorForReaction()Description: Retrieve the correction factor applied to a reaction in a zone. Syntax:
double
Libnucnet__Zone__getReverseRatioCorrectionFactorForReaction(
Libnucnet__Zone *self,
Libnucnet__Reaction *p_reaction
);
Input:
For valid input, the routine returns the correction factor that was applied in computing the rates for a given zone. If any input is invalid, error handling is invoked. Example: Print the correction factor that was applied to the reaction si28 + he4 -> s32 + gamma in the Libnucnet__Zone p_zone:
printf(
"Correction factor = %e\n",
Libnucnet__Zone__getReverseRatioCorrectionFactorForReaction(
p_zone,
Libnucnet__Reac__getReactionByString(
Libnucnet__Net__getReac(
Libnucnet__Zone__getNet( p_zone )
),
"si28 + he4 -> s32 + gamma"
)
)
);
|
|
|
Name: Libnucnet__Zone__getScreeningFactorForReaction()Description: Retrieve the screening factor applied to a reaction in a zone. Syntax:
double
Libnucnet__Zone__getScreeningFactorForReaction(
Libnucnet__Zone *self,
Libnucnet__Reaction *p_reaction
);
Input:
For valid input, the routine returns the screening factor that was applied in computing the rates for a given zone. If any input is invalid, error handling is invoked. Example: Print the screening factor that was applied to the reaction c12 + he4 -> o16 + gamma in the Libnucnet__Zone p_zone:
printf(
"Screening factor = %e\n",
Libnucnet__Zone__getScreeningFactorForReaction(
p_zone,
Libnucnet__Reac__getReactionByString(
Libnucnet__Net__getReac(
Libnucnet__Zone__getNet( p_zone )
),
"c12 + he4 -> o16 + gamma"
)
)
);
|
|
|
Name: Libnucnet__Zone__getSpeciesAbundance()Description: Retrieves the specified species abundance. Syntax:
double
Libnucnet__Zone__getSpeciesAbundance(
Libnucnet__Zone *self,
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__getNuc( p_my_nucnet ), "he4"
)
);
|
|
|
Name: Libnucnet__Zone__getSpeciesAbundanceChange()Description: Retrieves the specified species abundance change. Syntax:
double
Libnucnet__Zone__getSpeciesAbundanceChange(
Libnucnet__Zone *self,
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__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_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__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__Net__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__Net__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__Net__screeningFunction) my_screening_function,
&my_data
);
|
|
|
Name: Libnucnet__Zone__updateAbundanceChanges()Description: Routine to update the abundance changes of a given zone from an input gsl_vector. Syntax:
gsl_vector *
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:
gsl_vector *
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__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__updateAbundance( 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__updateAbundanceChange( p_zone, p_species, d_dy );
|
|
|
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:
Libnucnet__Zone *
Libnucnet__addZone(
Libnucnet *self,
int i_dim,
char *s_label1,
char *s_label2,
char *s_label3
);
Input:
Routine returns a pointer to the new zone which has also been added to the Libnucnet structure. If the zone could not be added, routine returns NULL. If any input is invalid, Libnucnet error handling is invoked. Example: Add a 1 dimensional zone with label "Z" to Libnucnet *p_my_nucnet:
Libnucnet__addZone( p_my_nucnet, 1, "Z" );
|
|
|
Name: Libnucnet__assignZonesFromXml()Description: Reads in initial mass fraction data from an XML file, creates zones, and stores the data in a Libnucnet structure. Syntax:
void
Libnucnet__assignZonesFromXml(
Libnucnet *self,
const char *s_xml_filename
);
Input:
If the input is valid, the Libnucnet structure has been updated with the zones. Example: Store the initial mass fraction data in input_mass.xml in p_my_nucnet:
Libnucnet__assignZonesFromXml( p_my_nucnet, "input_mass.xml" );
|
|
|
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( 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(
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__iterateZones()Description: Iterate through the zones and apply the user-supplied iterate function. Syntax:
void
Libnucnet__iterateZones(
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 in 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_xpath
);
Input:
If the input is valid, routine returns a pointer to a Libnucnet structure containing reactions selected by the xpath expression. If no reactions were found, p_reactions_list is empty. If the xpath expression 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 );
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", "[reactant = 'h1']"
);
|
|
|
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" )
);
|
|
|