Cantera  3.1.0b1
Loading...
Searching...
No Matches
Inlet1D Class Reference

An inlet. More...

#include <Boundary1D.h>

Inheritance diagram for Inlet1D:
[legend]

Detailed Description

An inlet.

Unstrained flows use an inlet to the left of the flow domain (left-to-right flow). Strained flow configurations may have inlets on the either side of the flow domain.

Definition at line 140 of file Boundary1D.h.

Public Member Functions

 Inlet1D ()
 Default constructor.
 
 Inlet1D (shared_ptr< Solution > solution, const string &id="")
 Constructor with contents.
 
string domainType () const override
 Domain type flag.
 
void setSpreadRate (double V0) override
 set spreading rate
 
double spreadRate () override
 Tangential velocity gradient [1/s] at this boundary.
 
void setTemperature (double T) override
 Set the temperature.
 
void show (const double *x) override
 Print the solution.
 
size_t nSpecies () override
 Get the number of species.
 
void setMoleFractions (const string &xin) override
 Set the mole fractions by specifying a string.
 
void setMoleFractions (const double *xin) override
 Set the mole fractions by specifying an array.
 
double massFraction (size_t k) override
 Mass fraction of species k.
 
void init () override
 Initialize.
 
void eval (size_t jg, double *xg, double *rg, integer *diagg, double rdt) override
 Evaluate the residual function at point j.
 
shared_ptr< SolutionArrayasArray (const double *soln) const override
 Save the state of this domain as a SolutionArray.
 
void fromArray (SolutionArray &arr, double *soln) override
 Restore the solution for this domain from a SolutionArray.
 
- Public Member Functions inherited from Boundary1D
 Boundary1D ()
 Default constructor.
 
void init () override
 Initialize.
 
string domainType () const override
 Domain type flag.
 
bool isConnector () override
 True if the domain is a connector domain.
 
virtual void setTemperature (double t)
 Set the temperature.
 
virtual double temperature ()
 Temperature [K].
 
virtual size_t nSpecies ()
 Get the number of species.
 
virtual void setMoleFractions (const string &xin)
 Set the mole fractions by specifying a string.
 
virtual void setMoleFractions (const double *xin)
 Set the mole fractions by specifying an array.
 
virtual double massFraction (size_t k)
 Mass fraction of species k.
 
virtual void setMdot (double mdot)
 Set the total mass flow rate [kg/m²/s].
 
virtual void setSpreadRate (double V0)
 Set tangential velocity gradient [1/s] at this boundary.
 
virtual double spreadRate ()
 Tangential velocity gradient [1/s] at this boundary.
 
virtual double mdot ()
 The total mass flow rate [kg/m2/s].
 
void setupGrid (size_t n, const double *z) override
 called to set up initial grid, and after grid refinement
 
void fromArray (SolutionArray &arr, double *soln) override
 Restore the solution for this domain from a SolutionArray.
 
- Public Member Functions inherited from Domain1D
 Domain1D (size_t nv=1, size_t points=1, double time=0.0)
 Constructor.
 
 Domain1D (const Domain1D &)=delete
 
Domain1Doperator= (const Domain1D &)=delete
 
virtual string domainType () const
 Domain type flag.
 
string type () const
 String indicating the domain implemented.
 
size_t domainIndex ()
 The left-to-right location of this domain.
 
virtual bool isConnector ()
 True if the domain is a connector domain.
 
void setSolution (shared_ptr< Solution > sol)
 Set the solution manager.
 
virtual void setKinetics (shared_ptr< Kinetics > kin)
 Set the kinetics manager.
 
virtual void setTransport (shared_ptr< Transport > trans)
 Set transport model to existing instance.
 
const OneDimcontainer () const
 The container holding this domain.
 
void setContainer (OneDim *c, size_t index)
 Specify the container object for this domain, and the position of this domain in the list.
 
void setBandwidth (int bw=-1)
 Set the Jacobian bandwidth. See the discussion of method bandwidth().
 
size_t bandwidth ()
 Set the Jacobian bandwidth for this domain.
 
virtual void init ()
 Initialize.
 
virtual void setInitialState (double *xlocal=0)
 
virtual void setState (size_t point, const double *state, double *x)
 
virtual void resetBadValues (double *xg)
 When called, this function should reset "bad" values in the state vector such as negative species concentrations.
 
virtual void resize (size_t nv, size_t np)
 Resize the domain to have nv components and np grid points.
 
Refinerrefiner ()
 Return a reference to the grid refiner.
 
size_t nComponents () const
 Number of components at each grid point.
 
void checkComponentIndex (size_t n) const
 Check that the specified component index is in range.
 
void checkComponentArraySize (size_t nn) const
 Check that an array size is at least nComponents().
 
size_t nPoints () const
 Number of grid points in this domain.
 
void checkPointIndex (size_t n) const
 Check that the specified point index is in range.
 
void checkPointArraySize (size_t nn) const
 Check that an array size is at least nPoints().
 
virtual string componentName (size_t n) const
 Name of component n. May be overloaded.
 
void setComponentName (size_t n, const string &name)
 Set the name of the component n to name.
 
virtual size_t componentIndex (const string &name) const
 index of component with name name.
 
void setBounds (size_t n, double lower, double upper)
 Set the upper and lower bounds for a solution component, n.
 
void setTransientTolerances (double rtol, double atol, size_t n=npos)
 Set tolerances for time-stepping mode.
 
void setSteadyTolerances (double rtol, double atol, size_t n=npos)
 Set tolerances for steady-state mode.
 
double rtol (size_t n)
 Relative tolerance of the nth component.
 
double atol (size_t n)
 Absolute tolerance of the nth component.
 
double steady_rtol (size_t n)
 Steady relative tolerance of the nth component.
 
double steady_atol (size_t n)
 Steady absolute tolerance of the nth component.
 
double transient_rtol (size_t n)
 Transient relative tolerance of the nth component.
 
double transient_atol (size_t n)
 Transient absolute tolerance of the nth component.
 
double upperBound (size_t n) const
 Upper bound on the nth component.
 
double lowerBound (size_t n) const
 Lower bound on the nth component.
 
void initTimeInteg (double dt, const double *x0)
 Performs the setup required before starting a time-stepping solution.
 
void setSteadyMode ()
 Set the internally-stored reciprocal of the time step to 0.0, which is used to indicate that the problem is in steady-state mode.
 
bool steady ()
 True if in steady-state mode.
 
bool transient ()
 True if not in steady-state mode.
 
void needJacUpdate ()
 Set this if something has changed in the governing equations (for example, the value of a constant has been changed, so that the last-computed Jacobian is no longer valid.
 
virtual void eval (size_t j, double *x, double *r, integer *mask, double rdt=0.0)
 Evaluate the residual function at point j.
 
size_t index (size_t n, size_t j) const
 Returns the index of the solution vector, which corresponds to component n at grid point j.
 
double value (const double *x, size_t n, size_t j) const
 Returns the value of solution component n at grid point j of the solution vector x.
 
virtual void setJac (MultiJac *jac)
 
virtual shared_ptr< SolutionArrayasArray (const double *soln) const
 Save the state of this domain as a SolutionArray.
 
shared_ptr< SolutionArraytoArray (bool normalize=false) const
 Save the state of this domain to a SolutionArray.
 
virtual void fromArray (SolutionArray &arr, double *soln)
 Restore the solution for this domain from a SolutionArray.
 
void fromArray (const shared_ptr< SolutionArray > &arr)
 Restore the solution for this domain from a SolutionArray.
 
shared_ptr< Solutionsolution () const
 Return thermo/kinetics/transport manager used in the domain.
 
size_t size () const
 Return the size of the solution vector (the product of m_nv and m_points).
 
void locate ()
 Find the index of the first grid point in this domain, and the start of its variables in the global solution vector.
 
virtual size_t loc (size_t j=0) const
 Location of the start of the local solution vector in the global solution vector.
 
size_t firstPoint () const
 The index of the first (that is, left-most) grid point belonging to this domain.
 
size_t lastPoint () const
 The index of the last (that is, right-most) grid point belonging to this domain.
 
void linkLeft (Domain1D *left)
 Set the left neighbor to domain 'left.
 
void linkRight (Domain1D *right)
 Set the right neighbor to domain 'right.'.
 
void append (Domain1D *right)
 Append domain 'right' to this one, and update all links.
 
Domain1Dleft () const
 Return a pointer to the left neighbor.
 
Domain1Dright () const
 Return a pointer to the right neighbor.
 
double prevSoln (size_t n, size_t j) const
 Value of component n at point j in the previous solution.
 
void setID (const string &s)
 Specify an identifying tag for this domain.
 
string id () const
 Returns the identifying tag for this domain.
 
virtual void show (std::ostream &s, const double *x)
 Print the solution.
 
virtual void show (const double *x)
 Print the solution.
 
double z (size_t jlocal) const
 Get the coordinate [m] of the point with local index jlocal
 
double zmin () const
 Get the coordinate [m] of the first (leftmost) grid point in this domain.
 
double zmax () const
 Get the coordinate [m] of the last (rightmost) grid point in this domain.
 
void setProfile (const string &name, double *values, double *soln)
 Set initial values for a component at each grid point.
 
vector< double > & grid ()
 Access the array of grid coordinates [m].
 
const vector< double > & grid () const
 Access the array of grid coordinates [m].
 
double grid (size_t point) const
 
virtual void setupGrid (size_t n, const double *z)
 called to set up initial grid, and after grid refinement
 
virtual void _getInitialSoln (double *x)
 Writes some or all initial solution values into the global solution array, beginning at the location pointed to by x.
 
virtual double initialValue (size_t n, size_t j)
 Initial value of solution component n at grid point j.
 
virtual void _finalize (const double *x)
 In some cases, a domain may need to set parameters that depend on the initial solution estimate.
 
void forceFullUpdate (bool update)
 In some cases, for computational efficiency some properties (such as transport coefficients) may not be updated during Jacobian evaluations.
 
void setData (shared_ptr< vector< double > > &data)
 Set shared data pointer.
 

Protected Attributes

int m_ilr
 A marker that indicates whether this is a left inlet or a right inlet.
 
double m_V0 = 0.0
 The spread rate of the inlet [1/s].
 
size_t m_nsp = 0
 Number of species in the adjacent flow domain.
 
vector< double > m_yin
 inlet mass fractions
 
string m_xstr
 inlet mass fractions. Parsing deferred to init()
 
Flow1Dm_flow = nullptr
 the adjacent flow domain
 
- Protected Attributes inherited from Boundary1D
Flow1Dm_flow_left = nullptr
 Flow domain to the left of this boundary.
 
Flow1Dm_flow_right = nullptr
 
size_t m_left_nv = 0
 Flow domain to the right of this boundary.
 
size_t m_right_nv = 0
 Number of state vector components in right flow domain.
 
size_t m_left_nsp = 0
 Number of species in left flow domain.
 
size_t m_right_nsp = 0
 Number of species in right flow domain.
 
ThermoPhasem_phase_left = nullptr
 Thermo object used by left flow domain.
 
ThermoPhasem_phase_right = nullptr
 Thermo object used by right flow domain.
 
double m_temp = 0.0
 Temperature of the boundary.
 
double m_mdot = 0.0
 Mass flow rate at the boundary.
 
- Protected Attributes inherited from Domain1D
shared_ptr< vector< double > > m_state
 data pointer shared from OneDim
 
double m_rdt = 0.0
 Reciprocal of the time step.
 
size_t m_nv = 0
 Number of solution components.
 
size_t m_points
 Number of grid points.
 
vector< double > m_slast
 Solution vector at the last time step.
 
vector< double > m_max
 Upper bounds on solution components.
 
vector< double > m_min
 Lower bounds on solution components.
 
vector< double > m_rtol_ss
 Relative tolerances for steady mode.
 
vector< double > m_rtol_ts
 Relative tolerances for transient mode.
 
vector< double > m_atol_ss
 Absolute tolerances for steady mode.
 
vector< double > m_atol_ts
 Absolute tolerances for transient mode.
 
vector< double > m_z
 1D spatial grid coordinates
 
OneDimm_container = nullptr
 Parent OneDim simulation containing this and adjacent domains.
 
size_t m_index
 Left-to-right location of this domain.
 
size_t m_iloc = 0
 Starting location within the solution vector for unknowns that correspond to this domain.
 
size_t m_jstart = 0
 Index of the first point in this domain in the global point list.
 
Domain1Dm_left = nullptr
 Pointer to the domain to the left.
 
Domain1Dm_right = nullptr
 Pointer to the domain to the right.
 
string m_id
 Identity tag for the domain.
 
unique_ptr< Refinerm_refiner
 Refiner object used for placing grid points.
 
vector< string > m_name
 Names of solution components.
 
int m_bw = -1
 See bandwidth()
 
bool m_force_full_update = false
 see forceFullUpdate()
 
shared_ptr< Solutionm_solution
 Composite thermo/kinetics/transport handler.
 

Additional Inherited Members

- Protected Member Functions inherited from Boundary1D
void _init (size_t n)
 Initialize member variables based on the adjacent domains.
 
- Protected Member Functions inherited from Domain1D
virtual AnyMap getMeta () const
 Retrieve meta data.
 
virtual void setMeta (const AnyMap &meta)
 Retrieve meta data.
 

Constructor & Destructor Documentation

◆ Inlet1D() [1/2]

Inlet1D ( )

Default constructor.

Definition at line 83 of file Boundary1D.cpp.

◆ Inlet1D() [2/2]

Inlet1D ( shared_ptr< Solution solution,
const string &  id = "" 
)

Constructor with contents.

Parameters
solutionSolution representing contents of adjacent flow domain
idName used to identify this domain

Definition at line 87 of file Boundary1D.cpp.

Member Function Documentation

◆ domainType()

string domainType ( ) const
inlineoverridevirtual

Domain type flag.

Since
Starting in Cantera 3.1, the return type is a string.

Reimplemented from Boundary1D.

Definition at line 151 of file Boundary1D.h.

◆ setSpreadRate()

void setSpreadRate ( double  V0)
overridevirtual

set spreading rate

Reimplemented from Boundary1D.

Definition at line 96 of file Boundary1D.cpp.

◆ spreadRate()

double spreadRate ( )
inlineoverridevirtual

Tangential velocity gradient [1/s] at this boundary.

Reimplemented from Boundary1D.

Definition at line 157 of file Boundary1D.h.

◆ setTemperature()

void setTemperature ( double  t)
overridevirtual

Set the temperature.

Reimplemented from Boundary1D.

Definition at line 102 of file Boundary1D.cpp.

◆ show()

void show ( const double *  x)
overridevirtual

Print the solution.

Parameters
xPointer to the local portion of the system state vector

Reimplemented from Domain1D.

Definition at line 111 of file Boundary1D.cpp.

◆ nSpecies()

size_t nSpecies ( )
inlineoverridevirtual

Get the number of species.

Reimplemented from Boundary1D.

Definition at line 165 of file Boundary1D.h.

◆ setMoleFractions() [1/2]

void setMoleFractions ( const string &  xin)
overridevirtual

Set the mole fractions by specifying a string.

Reimplemented from Boundary1D.

Definition at line 127 of file Boundary1D.cpp.

◆ setMoleFractions() [2/2]

void setMoleFractions ( const double *  xin)
overridevirtual

Set the mole fractions by specifying an array.

Reimplemented from Boundary1D.

Definition at line 137 of file Boundary1D.cpp.

◆ massFraction()

double massFraction ( size_t  k)
inlineoverridevirtual

Mass fraction of species k.

Reimplemented from Boundary1D.

Definition at line 171 of file Boundary1D.h.

◆ init()

void init ( )
overridevirtual

Initialize.

This method is called by OneDim::init() for each domain once at the beginning of a simulation. Base class method does nothing, but may be overloaded.

Reimplemented from Boundary1D.

Definition at line 146 of file Boundary1D.cpp.

◆ eval()

void eval ( size_t  j,
double *  x,
double *  r,
integer *  mask,
double  rdt 
)
overridevirtual

Evaluate the residual function at point j.

If j == npos, evaluate the residual function at all points.

This function must be implemented in classes derived from Domain1D.

Parameters
[in]jGrid point at which to update the residual
[in]xState vector
[out]rresidual vector
[out]maskBoolean mask indicating whether each solution component has a time derivative (1) or not (0).
[in]rdtReciprocal of the timestep (rdt=0 implies steady-state.)

Reimplemented from Domain1D.

Definition at line 179 of file Boundary1D.cpp.

◆ asArray()

shared_ptr< SolutionArray > asArray ( const double *  soln) const
overridevirtual

Save the state of this domain as a SolutionArray.

Parameters
solnlocal solution vector for this domain
Todo:
Despite the method's name, data are copied; the intent is to access data directly in future revisions, where a non-const version will be implemented.
Since
New in Cantera 3.0.

Reimplemented from Domain1D.

Definition at line 267 of file Boundary1D.cpp.

◆ fromArray()

void fromArray ( SolutionArray arr,
double *  soln 
)
overridevirtual

Restore the solution for this domain from a SolutionArray.

Parameters
[in]arrSolutionArray defining the state of this domain
[out]solnValue of the solution vector, local to this domain
Since
New in Cantera 3.0.

Reimplemented from Boundary1D.

Definition at line 284 of file Boundary1D.cpp.

Member Data Documentation

◆ m_ilr

int m_ilr
protected

A marker that indicates whether this is a left inlet or a right inlet.

Definition at line 181 of file Boundary1D.h.

◆ m_V0

double m_V0 = 0.0
protected

The spread rate of the inlet [1/s].

Definition at line 184 of file Boundary1D.h.

◆ m_nsp

size_t m_nsp = 0
protected

Number of species in the adjacent flow domain.

Definition at line 186 of file Boundary1D.h.

◆ m_yin

vector<double> m_yin
protected

inlet mass fractions

Definition at line 187 of file Boundary1D.h.

◆ m_xstr

string m_xstr
protected

inlet mass fractions. Parsing deferred to init()

Definition at line 188 of file Boundary1D.h.

◆ m_flow

Flow1D* m_flow = nullptr
protected

the adjacent flow domain

Definition at line 189 of file Boundary1D.h.


The documentation for this class was generated from the following files: