DDSL: Digital Dynamic Simulation Library for C/C++,based on Graph Theory

This software can be redistributed under GNU Lesser General Public License.
And every source codes can be obtained through GitHub
Copylight (c) 2020 Shigeo Kobayashi. All rights reserved.
[==>Japanese]
OverviewFlags(Attributes)API referencesReferences

DDSL(Digital Dynamic Simulation Library for C/C++,based on Graph Theory) is the library to assist numerical simulation for dynamics and statics.
DDSL implements Newton-method for solving non-linear algebraic eqution systems for static simulation and Euler,Backward Euler,and Runge-Kutta methods for solving ordinary differential equation systems for dynamic simulation.
The users,although,need not use these methos directly.
With the assistance of Graph-theory(about the theory,see References.), DDSL automaticaly determines computation sequence,constructs non-linear equation and ordinary differential equation systems if necessary,and solve them.

Overview

DDS_PROCESSORDDS_VARIABLE
To utilize DDSL,the user must define at least one processor(
DDS_PROCESSOR) and more than or equal to one variable(DDS_VARIABLE).
The processor(DDS_PROCESSOR) keeps every variables(DDS_VARIABLE) and solves algebraic and ordinary differential equations.

As is shown in the examples bellow,setting DDS_FLAG_REQUIRED flag is the most important for the user. DDSL does everything,including computation order construction, to compute the variables of DDS_FLAG_REQUIRED flag.

DDS_PROCESSOR and DDS_VARIABLE are pointers of the structures defined in ddsl.h.

Note: Include ddsl.h in the source code and link appropriate libraries(including DDSL library).
Note: Any name offered by this library has prefix "DDS_" or "Dds".

DDS_PROCESSOR

To define DDS_PROCESSOR:
    DDS_PROCESSOR p;
    int e = DdsCreateProcessor(&p,10);
    if(e!=0) printf("DdsCreateProcessor() returned an error\n");
DdsCreateProcessor(&p,10) allocates memories for DDS_PROCESSOR structure and sets it's pointer to p. The second argument 10 is the number of variables(explained later) the DdsProcessor() can keep which will be automatically extended if necessary. DdsCreateProcessor() returns 0 if DdsProcessor is normally allocated,or returns non-zero error code.

Note:DdsCreateProcessor() automatically allocates default variables,time and integration step,internally.

DDS_VARIABLE in the simple static example

The variable(DDS_VARIABLE) has value(double),flags/attributes(unsigned int),and may have right hand side variables(hereinafter referred to as RHSVs) and the function pointer to compute it's value(Note: The user must not change the variable's value and any other's within that function. DDSL does that instead of the user.).
Definition:
Upstream: In case like y = f(u,v,w). Variables u,v,and w are upstream variables of y.
Downstream: In case like x = f(y). Variable x is a downstream variable of y.
Variables must be registered to the processor(DDS_PROCESSOR) as:
    DDS_PROCESSOR p;
    DDS_VARIABLE y,x1,x2;
    DdsCreateProcessor(&p,10);
    DdsAddVariableV(p,&x1, "x1", DDS_FLAG_SET, 1.0, NULL,0); /* Register x1 to p. */
    DdsAddVariableV(p,&x2, "x2", DDS_FLAG_SET, 2.0, NULL,0);
    DdsAddVariableV(p,&y,  "y",  DDS_FLAG_REQUIRED, 0.0, CompY,2,x1,x2);
    DdsCompileGraph(p, 0); /* Check relations and determine computation order. */
    DdsComputeStatic(p);   /* Compute variable's value according to defined(by DdsCompileGraph()) order. */
    printf("Value: y=%lf x1=%lf x2=%lf\n",DdsGetValue(y),DdsGetValue(x1),DdsGetValue(x2));
    /* ==> Value: y=0.718282 x1=1.000000 x2=2.000000 */
    DdsDeleteProcessor(&p);
   double CompY(DDS_PROCESSOR p, DDS_VARIABLE y)
   {
	DDS_VARIABLE x1 = DdsGetRHSV(y, 0);
	DDS_VARIABLE x2 = DdsGetRHSV(y, 1);
	double xv1 = DdsGetValue(x1);
	double xv2 = DdsGetValue(x2);
	return exp(xv1)-xv2;
   }
Note:Functions used in the example may return non-zero value in case of error,that is ignored to make the example be simple.

DdsAddVariableV(p,&y, "y", DDS_FLAG_REQUIRED, 0.0, CompY,2,x1,x2)
allocates memories of DDS_VARIABLE,sets it's pointer to y,and register the newly allocated variable(y) to the processor p.
Also DdsAddVariableV() sets flag DDS_FLAG_REQUIRED,value 0.0,function pointer(CompY) to compute the value of y,and right hand side variables x1,x2(2 is the number of right hand side variables, x1,x2), to y.
Note: Note that right hand side variables x1 and x2 are already registered to p before y is,and have valid address registered to p at this moment.

The flag DDS_FLAG_REQUIRED is the specification that the value of the variable y is what the user wants to know.

DdsAddVariableV(p,&x1, "x1", DDS_FLAG_SET, 1.0, NULL,0)
The flag DDS_FLAG_SET for x1 and x2 is a constant requirement,and the value specified (1.0 or 2.0) will never be changed even right hand variables are given.

DdsCompileGraph(p, 0) checks relations and determines computation order,and DdsComputeStatic(p) simply computes the required value of y.

Instead of redundancy above,the user will obtain flexibility.
If the flags are changed to as:

    DdsAddVariableV(p,&x1, "x1", DDS_FLAG_REQUIRED, 1.0, NULL,0);
    DdsAddVariableV(p,&y,  "y",  DDS_FLAG_TARGETED, 0.0, CompY,2,x1,x2);
    /* ==> Value: y=0.000000 x1=0.693148 x2=2.000000 */
DDS_FLAG_TARGETED means that the value of y is to be 0.0 after computations. And the value of x1,which is not DDS_FLAG_SET this time,will be the value that make the value of y be 0.0. DDSL creates (non-linear) equation(s) and solves it automatically according to the user's specifications.

Note: Flags can be freely changed after DdsAddVariableV().

DDS_FLAG_INTEGRATED and DDS_STEADY_STATE dynamic example

The next example demonstrates the integration(dynamic) and the computation of the steady state.
	/* dy/dt = A(C-y) + B */
	/* dy/dt(inf)==>0.0, y = (AC+B)/A  */
	/* Steady state y = 5.0 */
	DDS_PROCESSOR p;
	DDS_VARIABLE y, dydt, A,B,C, time, step;
	DDS_VARIABLE* pVs, *pVr;
	int i,j,nv, nr;
	DdsCreateProcessor(&p, 10);
	DdsAddVariableV(p, &A, "A", DDS_FLAG_SET, 1.0, NULL, 0);
	DdsAddVariableV(p, &B, "B", DDS_FLAG_SET, 2.0, NULL, 0);
	DdsAddVariableV(p, &C, "C", DDS_FLAG_SET, 3.0, NULL, 0);
	DdsAddVariableV(p, &y, "y", DDS_FLAG_REQUIRED|DDS_FLAG_INTEGRATED, 1.0, NULL,1, &dydt);
	DdsAddVariableV(p, &dydt, "dy/dt", DDS_FLAG_REQUIRED, 0.0, CompDYDT, 4,&A,&B,&C,&y);
	time = DdsTime(p);
	step = DdsStep(p);
	pVs = DdsGetVariables(&nv, p);
	for (i = 0; i < nv; ++i) {
		pVr = DdsGetRhsvs(&nr, pVs[i]);
		for (j = 0; j < nr; ++j) {
			pVr[j] = *((DDS_VARIABLE*)pVr[j]);
		}
	}
	DdsSetValue(step, 0.1);
	DdsCompileGraph(p, 0);
	for (i = 0; i < 11; ++i) {
		DdsComputeDynamic(p,0);
		printf("Dynamic:Time=%lf dydt=%lf y=%lf\n",
			DdsGetValue(time), DdsGetValue(dydt), DdsGetValue(y));
	}
	DdsCompileGraph(p, DDS_STEADY_STATE);
	DdsComputeStatic(p);
	printf("Steady :Time=%lf dydt=%lf y=%lf\n",
		DdsGetValue(time), DdsGetValue(dydt), DdsGetValue(y));
	DdsDeleteProcessor(&p);
double CompDYDT(DDS_PROCESSOR p, DDS_VARIABLE y)
{
	DDS_VARIABLE A = DdsGetRHSV(y, 0);
	DDS_VARIABLE B = DdsGetRHSV(y, 1);
	DDS_VARIABLE C = DdsGetRHSV(y, 2);
	DDS_VARIABLE Y = DdsGetRHSV(y, 3);
	double AV = DdsGetValue(A);
	double BV = DdsGetValue(B);
	double CV = DdsGetValue(C);
	double yv = DdsGetValue(Y);
	return AV*(CV-yv)+BV;
}

/* Computation results 
Dynamic:Time=0.000000 dydt=4.000000 y=1.000000
Dynamic:Time=0.100000 dydt=3.619000 y=1.380650
Dynamic:Time=0.200000 dydt=3.274607 y=1.725071
Dynamic:Time=0.300000 dydt=2.962992 y=2.036717
Dynamic:Time=0.400000 dydt=2.681031 y=2.318706
Dynamic:Time=0.500000 dydt=2.425901 y=2.573860
Dynamic:Time=0.600000 dydt=2.195050 y=2.804734
Dynamic:Time=0.700000 dydt=1.986167 y=3.013638
Dynamic:Time=0.800000 dydt=1.797161 y=3.202662
Dynamic:Time=0.900000 dydt=1.626141 y=3.373699
Dynamic:Time=1.000000 dydt=1.471396 y=3.528459
Steady :Time=1.000000 dydt=0.000000 y=5.000000
*/
Any variable including RHSVs of any variable must be registered to the processor before calling DdsCompileGraph(). But this is sometimes difficult in programing(when DdsAddVariableV() is called),especially the case of loop.
The following codes(extracted from above) avoid the difficulties by specifying pointers for RHSVs at first,and reset them as follows(see API references):
	pVs = DdsGetVariables(&nv, p);
	for (i = 0; i < nv; ++i) {
		pVr = DdsGetRhsvs(&nr, pVs[i]);
		for (j = 0; j < nr; ++j) {
			pVr[j] = *((DDS_VARIABLE*)pVr[j]);
		}
	}
Note: All RHSVs of any variable are all pointers this time(ex:DdsAddVariableV(p, &dydt, "dy/dt", DDS_FLAG_REQUIRED, 0.0, CompDYDT, 4,&A,&B,&C,&y)).
Note: DdsAddVariableV() never check any RHSV given(so RHSV can even be NULL).

DdsAddVariableV(p, &y, "y", DDS_FLAG_REQUIRED|DDS_FLAG_INTEGRATED, 1.0, NULL,1, &dydt)
means that the y,which has DDS_FLAG_INTEGRATED flag,is computed by dydt through integration.
As is shown in the left image,all variables other than y are computed at the specific time(on the specific time plain).
y,on the other hand, is computed(integrated) when time is incremented.
DdsComputeDynamic(p,0) integrates all variables of DDS_FLAG_INTEGRATED by specified step(1 step), and after integration, it calls DdsComputeStatic(p) internally and increment the time by specified step.

Note: DDSL can compute the steady state by calling DdsCompileGraph(p, DDS_STEADY_STATE).
Note: y(=DDS_FLAG_INTEGRATED:<I>) looks like constant(=DDS_FLAG_SET:<S>) from downstream variables of y(see flags).
Note: More explanation will be given in detail hearinafter.

What the user must do:

As described so far,the user must invoke:
  1. DdsCreateProcessor()...create a processor
  2. DdsAddVariableV()
    DdsAddVariableA()...create a variable and register it to the processor. The user also sets the (initial) value of y,the right hand side variables,and the function to compute the value of y. The user also sets user FLAGs(attributes) of y that determin the computation order.
  3. DdsCompileGraph()... checks relations of each variable and determines computation order.
  4. DdsComputeStatic() or
    DdsComputeDynamic()... compute variable's value by solving non-linear equations or ordinary differential equations,according to the computation order determined by DdsCompileGraph().
  5. DdsDeleteProcessor()... frees every memories allocated by the processor.

What DDSL can do:

DDSL can:
  1. remove variables having no contribution directly or indirectly for computing variables of DDS_FLAG_REQUIRED from any computation.
  2. automatically detect loop such as y=f(y),and rearrange it to the equation system as y2 = y1 - f(y1), y2 ==>0.0
  3. automatically construct non-linear eqution sysytems and solve them.
  4. construct computation order of each variable.
  5. solve ordinally differential equations if given.
  6. obtain steady state by changing ordinally differential equations to algebraic equation systems and solve them.
How DDSL realizes these will be hereinafter described.

Flags(Attributes)

User flags System flags Block triangular decomposition
The flags or attributes (see ddsl.h) are the most impotant part. A variable has two 32-bit unsigned int element for user's flags and system's flags. User's flags are controled by the user and DDSL never change them. One flag occupies 1-bit of 32-bit unsigned integer which can be set or reset. DDSL checks user's flags and copies them to system's if they are correct. DDSL controls(set or reset) system's flag during the processing.
More than one flags can be specified as:

    DdsAddVariableV(p,&x2, "x2", DDS_FLAG_SET|DDS_FLAG_VOLATILE, 2.0, NULL,0);
or

    DdsSetUserFlag(x2,DDS_FLAG_SET|DDS_FLAG_VOLATILE);

User flags:

DDS_FLAG_REQUIREDDDS_FLAG_SETDDS_FLAG_TARGETEDDDS_FLAG_INTEGRATED
DDS_FLAG_VOLATILEDDS_FLAG_DIVISIBLEDDS_FLAG_NON_DIVISIBLE
Note: All user flags use lower 12-bit of unsigned int,thus,upper 20-bit can be freely used by the user.

User flag name (notation)value
DDS_FLAG_REQUIRED (<R>)0x00000001
The user must set this flag to any variable of which value is wanted.
If there is no variable of <R> flag,then nothing is computed.

DDS_FLAG_SET (<S>)0x00000002
<S> is a constant specification.
Any variable having this flag is never computed even RHSVs and the function to compute its value are defined.
This specification may cause to divide one connected component into more than one.
Note:Connected component is a group of variables that are connected directly or indirectly each other.

DDS_FLAG_TARGETED (<T>)0x00000004
The value of <T> variable becomes specified value by computation.
The variable which can freely change it's value(<F>(detected by DDSL) variable which has no RHSVs and not <S>) must exist on the upstream of <T>. The value of <F> variable will be adjusted to make the value of <T> be the specified value. The (computed) value of <T> variable becomes the specified value by solving non-linear equation: y<T> = f(x<F>).

Solvability conditions on non-linear equation systems(checked by DDSL):

  1. Number of <T>s and <F>s must be equal in the same connected component.
  2. Independen route from each <F> to <T> must exist(1 to 1 matching of <F> and <T>).
Note: The value of <F> before computation will be used as an initial value.
Note: N-<T>s(N-<F>s) construct N-dimensional equation system in general(can be reduced if possible) if they are mutually connected.

From downstream variables,<T> look like <S>.
This specification may cause to divide one connected component into more than one.


This looks to be 2-dimensional equation system because there are two <F>s and two <T>s. But,the independent 2 routes from each <F> to <T> can not be found that means the equation system can not be solved which is detected before actual computations(DDSL checks this condition).

DDS_FLAG_DIVISIBLE (<D>)0x00000020
Any variable(say y) computed by itself like y = f(y) constructs a loop. It is necessary to divide y ( y ==> y1 , y2 ) to compute value of the variables on the loop. If y is divided into y1(=<F>) and y2(=<T>),then a new equation
y2 = y1 - f(y1)
y2==>0.0
will be added and be solved.
y1 is <F> and y2 is <T>.also both y1 and y2 are <DD> at the same time.

DDSL automatically detects loops and performs above procedures.
Any variable on a loop can be divided. The user can specify this <D>(DDS_FLAG_DIVISIBLE) flag on any variable(by reason of easy estimation of initial value, etc.) to be divided. DDSL tries to divide <D> varable,but it is not alway divided. DDSL finds loops and tries to divide the variable which is located on more loops. The variable actually divided will have the system flag DDS_SFLAG_DIVIDED(<DD>).

Note: <D> variable is not alway divided.
Note: The y1 is the same as y,and newly created variable y2 has the name of y+'+'.

DDS_FLAG_NON_DIVISIBLE (<N>)0x00000040
This flag is set on the variable which the user do not want to be divided. DDSL tries not to divide the variable of this flag,but it is not alway true due to
above reason.

DDS_FLAG_INTEGRATED (<I>)0x00000008
Integration specification. <I> variable's value is computed by integration. <I> must have only 1 RHSV which is the derivative(
<DR>)of the ordinary differential equation. The function pointer to compute(integrate) the value of <I> is not necessary because integration methods such as Euler,Backward Euler,and Runge-Kutta methods are built in.

Steady state computation:

The steady state,if exists,can be obtained to keep integration which is,althogh,time consuming.

DDSL can algebraicaly compute the steady state as:

  1. Change <I> variables to <F>.
  2. Change <DR> variables to <T> of value 0.0(steady state condition).
  3. Solve non-linear equation systems to find the value of <I>s that make all value of <DR>s be zero(steady state).
See: DdsCompileGraph(p,DDS_STEADY_STATE)
DDS_FLAG_VOLATILE (<V>)0x00000010
<V> flag can only be set on the variable of <S> or <T>.
<V> flag means that the value of <V> is outside the control of DDSL(<V> treated as <S> or <T>,but it's value is undefined).
<V> variables are regarded as time-dependent variables.

Block triangular decomposition

Newton method creates Jacobian matrix when it solves non-linear equation system (y=f(x) ==>0.0),
y(x+dx) = y(x)+J*dx ==> 0.0 ,
then dx = -J-1*y(x),
where J is a Jacobian matrix(=dy/dx), y,x,and dx are vectors,and x+dx is the next estimation.
Computation of J or J-1 may be time-consuming if J is quite large.
DDSL divides,if possible,J into a sequence of small-sized matrices as shown in the left image which is called as Block Triangular Decomposition.
In the left image,instead of solving big J,smaller sized matrices (1),(2),and (3) are solved in this order by re-arranging columns and rows,
DDSL adopts numerical differentiation,so that the user need not to prepaire any equations to compute the Jacobian-matrix(dy/dx).

Systems flags:

DDS_SFLAG_ALIVEDDS_SFLAG_FREEDDS_SFLAG_DIVIDEDDDS_SFLAG_DERIVATIVE
DDS_SFLAG_ERRORDDS_COMPUTED_ONCEDDS_COMPUTED_EVERY_TIMEDDS_COMPUTED_ANY_TIME

As descrived before,DDSL copies user flags to system flag area after checking,and adds or changes them during the processing(DDSL never touch user flag area which is unsigned int).Typical system flags are listed and described here.

System flag name (notation)value
DDS_SFLAG_ALIVE (<AL>)0x00001000
Any variable which has,directly or indirectly,any contribution to compute <R> variables will be <AL>(variable aliving).
On the other hand,variables not <AL> will be removed from any computation.

DDS_SFLAG_FREE (<F>)0x00002000
DDSL finds and sets <F>flag on the variable which has no RHSV and not <S>.
<F> variables will be independent(because it can freely change it's value) variables of non-linear eqution systems.
Corresponding dependent variables are
<T>s(that must have <F>s on their upstream).

Note: Total number of <F>s and <T>s must be equal in the same connected component.
Note: Independent route from each <F> to <T> must exist.

DDS_SFLAG_DIVIDED (<DD>)0x00004000
DDSL detects a loop and divide a variable on the loop.
Dividing a variable means that one more variable is created(and registered to the processor).
The variable before divided is changed to <F>(and is <DD> at the same time),and the newly created variable is <T>(and also is <DD> at the same time).

DDS_SFLAG_DERIVATIVE (<DR>)0x00008000
<DR> is the RHSV of <I> and is the derivative of the <I>.
On steady state computation,<DR> will be changed to <T> of the value 0.0 (<I> will be changed to <F> at the same time).

DDS_SFLAG_ERROR (<ER>)0x00010000
When combination of user flags have any error,then this flag will be set to the system flag.

DDS_COMPUTED_ONCE (<1T>)0x00020000
Any variable which is computed only by non-
volatile <S> or non-volatile <T> directly or indirectly. Once <1T> is computed then it's value is unchanged.
Note: See DdsGetVariableSequence()

DDS_COMPUTED_EVERY_TIME (<ET>)0x00040000
<ET> variable is computed at every integration step. <ET> variable is also on the route from <I> to <DR>.
Note:
See DdsGetVariableSequence()

DDS_COMPUTED_ANY_TIME (<AT>)0x00080000
<AT> variable is time dependent but is not necessary to be computed at every integration step. If the user wants to know the value of <AT> variable,user must call
DdsComputeStatic() after DdsComputeDynamic().
Note: See DdsGetVariableSequence()

API references

Basic Processor FunctionsBasic Variable FunctionsService FunctionsLow Level Functions
DdsCreateProcessor
DdsDeleteProcessor
DdsAddVariableV
DdsAddVariableA
DdsCompileGraph
DdsComputeStatic
DdsComputeDynamic
DdsGetVariables
DdsTime
DdsStep
DdsGetVariableSequence
DdsSetErrorHandler
DdsGetErrorHandler
DdsGetVariableName
DdsGetValue
DdsSetValue
DdsGetRhsvs
DdsSetRHSV
DdsGetRHSV
DdsGetUserFlag
DdsSetUserFlag
DdsGetSystemFlag
DdsSetUserFlagOn
DdsSetUserFlagOff
DdsGetEPS
DdsSetEPS
DdsGetMaxIterations
DdsSetMaxIterations
DdsGetProcessorUserPTR
DdsSetProcessorUserPTR
DdsGetVariableUserPTR
DdsSetVariableUserPTR
DdsCheckVariable
DdsDbgPrintF
DdsGetVariableNext
DdsGetVariableIndex
DdsGetVariableScore
DdsSieveVariable
DdsDivideLoop
DdsCheckRouteFT
DdsBuildSequence

Return valueFunction name(arguments...)
intDdsCreateProcessor(DDS_PROCESSOR* p, int nv)
DdsCreateProcessor() allocates memories of DDS_PROCESSOR (struct defined ddsl.h) and sets it's pointer to the first argument p.
The second argument nv is the number of variables the processor p can keep which will be extended automaticaly if more than nv variables are to be added.
Note: This function returns 0 if it normally finishes,otherwise it returns non-zero error value(see error code listed in ddsl.h).

voidDdsDeleteProcessor(DDS_PROCESSOR* p)
DdsDeleteProcessor() frees all memories the processor p allocated.
And sets p be NULL.

intDdsAddVariableV(DDS_PROCESSOR p,DDS_VARIABLE *pv,const char *name,unsigned int f,double val, ComputeVal func,int nr,...)
DdsAddVariableV() allocates memories of DDS_VARIABLE (struct defined in ddsl.h) and sets it's pointer to pv.
And register newly created variable to the processor p. [Arguments]

Note: The last nr RHSVs of pv's can be anything(can be NULL) at this moment. DDSL allocates the memories for them and just stores them(see).
Note: This function returns 0 if it normally finishes,otherwise it returns non-zero error value(see error code listed in ddsl.h).

intDdsAddVariableA(DDS_PROCESSOR p, DDS_VARIABLE* pv, const char* name, unsigned int f, double val, ComputeVal func, int nr, DDS_VARIABLE* rhsvs)
DdsAddVariableA() is almost the same as
DdsAddVariableV() above except for the last argument.
The last argument of DdsAddVariableA() is an array of RHSVs of size nr.
Note: This function returns 0 if it normally finishes,otherwise it returns non-zero error value(see error code listed in ddsl.h).

intDdsCompileGraph(DDS_PROCESSOR p,int method)
DdsCompileGraph() calls following functions in this order and build up everything necessary before computation.

  1. DdsSieveVariable(p)... copies correct user flags to system flag,and selects <AL> variables.
  2. DdsDivideLoop(p)... detects loops and divide them into <DD>.
  3. DdsCheckRouteFT(p)...checks 1 to 1 correspondense between <F>s and <T>s,and performs block triangular decomposition.
  4. DdsBuildSequence(p)... constructs conputation order.
DdsCompileGraph() also divide all variables computed into 3 groups(DdsBuildSequence(p)) as: The method must be 0(default:DDS_I_RUNGE_KUTTA) or the one listed bellow: The computation order changes according to integration methods or steady state computation.

Note: This function returns 0 if it normally finishes,otherwise it returns non-zero error value(see error code listed in ddsl.h).

intDdsComputeStatic(DDS_PROCESSOR p)
DdsComputeStatic() computes value of variables at specific time,solves non-linear equation systems if given.
Note: This function returns 0 if it normally finishes,otherwise it returns non-zero error value(see error code listed in ddsl.h).

intDdsComputeDynamic(DDS_PROCESSOR p,int method)
DdsComputeDynamic() solves(1 step of) ordinaly differential equations by specified integration method(Euler,Backword Euler,or Runge-Kutta method).
DdsComputeDynamic() performs processes in the following order:

  1. DdsComputeDynamic() calls DdsComputeStatic() at t(=current time),if necessary.
  2. Suppose y is an <I>,then the value y(t+dt) is computed by integration(t is the current time and dt is the time step).
  3. DdsComputeDynamic(),then,calls DdsComputeStatic() to compute the value of any other time dependent variables(see: DdsGetVariableSequence(p,DDS_COMPUTED_EVERY_TIME) at time t+dt,and increments time = t + dt.
method can be 0(default), DDS_I_EULER, or DDS_I_RUNGE_KUTTA.
DDS_I_EULER or DDS_I_RUNGE_KUTTA can be specified only when the second argument of DdsCompileGraph() is the one of DDS_I_EULER or DDS_I_RUNGE_KUTTA(mutualy exchangeable).
Note: To increment Time automaticaly,the value of Step must be defined beforehand.
Note: This function returns 0 if it normally finishes,otherwise it returns non-zero error value(see error code listed in ddsl.h).

DDS_VARIABLEDdsTime(DDS_PROCESSOR p)
DdsTime() returns the builtin variable of which name is "#TIME" representing the global time.
The initial value is 0.0,and is automaticaly incremented by
DdsComputeDynamic()

DDS_VARIABLEDdsStep(DDS_PROCESSOR p)
DdsStep() returns the builtin variable of which name is "#STEP" representing the global time step.
The initial value is 0.0. User must define proper value of this before calling
DdsComputeDynamic()(see DdsSetValue()).

unsigned intDdsGetUserFlag(DDS_VARIABLE v)
DdsGetUserFlag() returns
every user flags.

unsigned intDdsSetUserFlag(DDS_VARIABLE v,unsigned int f)
DdsSetUserFlag() replaces
every user flags by f and returns f it self.

unsigned intDdsGetSystemFlag(DDS_VARIABLE v)
DdsGetSystemFlag() returns
every system flags.

unsigned intDdsSetUserFlagOn(DDS_VARIABLE v, unsigned int f)
DdsSetUserFlagOn() sets user flags specified by f be on as.
F |= f
where F is the user flags.

unsigned intDdsSetUserFlagOff(DDS_VARIABLE v, unsigned int f)
DdsSetUserFlagOff() sets user flags specified by f be off as.
F &= ~(f)
where F is the user flags.

DDS_VARIABLE*DdsGetVariables(int* nv, DDS_PROCESSOR p)
DdsGetVariables() returns the array of all variables registered to the processor p.
nv is the array size.

void*DdsGetVariableUserPTR(DDS_VARIABLE v)
DdsGetVariableUserPTR() returns the pointer(void*) set to v by the user. The usage of the pointer is decided by the user,DDS never touch this.

voidDdsSetVariableUserPTR(DDS_VARIABLE v, void* val)
DdsSetVariableUserPTR() sets the pointer(void*) val to to v. The usage of the pointer is decided by the user,DDS never touch this.

const char*DdsGetVariableName(DDS_VARIABLE v)
DdsGetVariableName() returns the name of the variable v.

doubleDdsGetValue(DDS_VARIABLE v)
DdsGetValue() returns the value of y.

doubleDdsSetValue(DDS_VARIABLE v,double val)
DdsSetValue() sets the value of y to val.

DDS_VARIABLE*DdsGetRhsvs(int* nr, DDS_VARIABLE v)
DdsGetRhsvs() returns the array of all right hand side variables registered to the variable v. nr is the number of the right hand side variables of v(array size).

intDdsSetRHSV(DDS_VARIABLE v, int i, DDS_VARIABLE rv)
DdsSetRHSV() sets i-th right hand side variable of v be rv.

DDS_VARIABLEDdsGetRHSV(DDS_VARIABLE v, int i)
DdsGetRHSV() returns i-th right hand side variable of v.

DDS_VARIABLEDdsGetVariableSequence(DDS_PROCESSOR p, unsigned int seq)
DDS determines computation order of each variable(see
DdsCompileGraph()).
DDS(DdsCompileGraph()) also make group of variables in which variables are linked according to the computation order as:

  1. DDS_COMPUTED_ONCE: Variable in this group is computed at once at the beginning and remain unchanged.
  2. DDS_COMPUTED_EVERY_TIME: Variable in this group is computed at every time step.
  3. DDS_COMPUTED_ANY_TIME: Variable in this group is computed(by DdsComputeStatic()) at any time after DdsComputeDynamic().
According to the computation order,variables in above group can be listed like as follows.
The argument seq must be the one listed above.:

  DDS_VARIABLE v = DdsGetVariableSequence(p,DDS_COMPUTED_ONCE);
  /* The order of variables computed at DDS_COMPUTED_ONCE */
  while(v != NULL) {
     printf(" Name=%s\n",DdsGetVariableName(v));
     v = DdsGetVariableNext(v);
  }
ErrHandlerDdsGetErrorHandler(DDS_PROCESSOR p)
DdsGetErrorHandler() returns ErrorHandler if registered to the processor p.
The user can obtain some error informations through this function on error(see ddsl.h).

voidDdsSetErrorHandler(DDS_PROCESSOR p, ErrHandler handler)
DdsSetErrorHandler() sets ErrorHandler(handler) to the processor p.
The user can obtain some error informations through this function on error(see ddsl.h).

doubleDdsGetEPS(DDS_PROCESSOR p)
DdsGetEPS() returns convergence criteria for solving non-linear equation systems.

doubleDdsSetEPS(DDS_PROCESSOR p, double eps)
DdsSetEPS() sets convergence criteria to eps for solving non-linear equation systems.

intDdsGetMaxIterations(DDS_PROCESSOR p)
DdsGetMaxIterations() returns the maximum iteration count for solving non-linear equation systems.

intDdsSetMaxIterations(DDS_PROCESSOR p, int max)
DdsSetMaxIterations() sets the maximum iteration count to max for solving non-linear equation systems.

void*DdsGetProcessorUserPTR(DDS_PROCESSOR p)
DdsGetProcessorUserPTR() returns the pointer(void*) set by the user. The usage of the pointer is decided by the user,DDS never touch this.

voidDdsSetProcessorUserPTR(DDS_PROCESSOR p,void* val)
DdsSetProcessorUserPTR() sets the pointer(void*) val to to p. The usage of the pointer is decided by the user,DDS never touch this.

intDdsCheckVariable(DDS_PROCESSOR p, DDS_VARIABLE hv)
DdsCheckVariable() checks the validity of the variables registered to p,varidity of user flags for instance.
If any error is found,then
DDS_SFLAG_ERROR is set to the variable of error.
DdsCheckVariable() returns non-zero value when any error is detected.

voidDdsDbgPrintF(FILE* f, const char* title, DDS_PROCESSOR p)
DdsDbgPrintF() outputs current variable informations to f.
For details,see the souce code ddsl.cpp.

DDS_VARIABLEDdsGetVariableNext(DDS_VARIABLE v)
DdsGetVariableNext() returns 'NEXT' variable.
Each variables has a link to the next variable of which usage depends on the process context.

intDdsGetVariableIndex(DDS_VARIABLE v)
DdsGetVariableIndex() returns 'INDEX' which is an integer.
Each variables has a integer value of which usage depends on the process context.

intDdsGetVariableScore(DDS_VARIABLE v)
Each variables has one more integer value of which usage also depends on the process context.

intDdsSieveVariable(DDS_PROCESSOR p)
See
DdsCompileGraph() source code.

intDdsDivideLoop(DDS_PROCESSOR p)
See
DdsCompileGraph() source code.

intDdsCheckRouteFT(DDS_PROCESSOR p)
See
DdsCompileGraph() source code.

intDdsBuildSequence(DDS_PROCESSOR p)
See
DdsCompileGraph() source code.


References:

Shigeo Kobayashi 2021-6-10

I dedicate this software to two great scientists&mathematicians,my former boss Junkichi Tunekawa and the professor Masao Iri. (They are no longer with us. This software is the proof that I had a time,although it was very short period,to work with them.)