13 #ifndef __deal2__appframe__DAE_solver_h
14 #define __deal2__appframe__DAE_solver_h
22 #define omp_get_thread_num() 0
27 void coldae_(
int &,
int &,
int [],
double &,
double &,
28 double [],
int [],
int [],
double [],
double [],
int [],
30 void (*)(
double &,
double [],
double [],
double []),
31 void (*)(
double &,
double [],
double [],
double []),
32 void (*)(
int &,
double [],
double &),
33 void (*)(
int &,
double [],
double []),
34 void (*)(
double &,
double [],
double [],
double []));
35 void appsln_ (
double &,
double [],
double [],
double [],
int []);
42 namespace ApplicationCore
46 typedef void (*
fsub_ptr)(
double &,
double [],
double [],
double []);
47 typedef void (*
dfsub_ptr)(
double &,
double [],
double [],
double []);
48 typedef void (*
gsub_ptr)(
int &,
double [],
double &);
49 typedef void (*
dgsub_ptr)(
int &,
double [],
double []);
50 typedef void (*
guess_ptr)(
double &,
double [],
double [],
double []);
124 DAESolver(
int m_comp,
int ny,
int m[],
double a,
double b,
double zeta[],
125 void (*fsub)(
double &,
double [],
double [],
double []),
126 void (*dfsub)(
double &,
double [],
double [],
double []),
127 void (*gsub)(
int &,
double [],
double &),
128 void (*dgsub)(
int &,
double [],
double []),
129 void (*guess)(
double &,
double [],
double [],
double [])=NULL);
265 void set_fsub(
void (*fsub)(
double &,
double [],
double [],
double []));
270 void set_dfsub(
void (*dfsub)(
double &,
double [],
double [],
double []));
275 void set_gsub(
void (*)(
int &,
double [],
double &));
280 void set_dgsub(
void (*)(
int &,
double [],
double []));
285 void set_guess(
void (*)(
double &,
double [],
double [],
double []));
289 void (*
DAE_fsub) (
double &,
double [],
double [],
double []);
291 void (*
DAE_dfsub) (
double &,
double [],
double [],
double []);
297 void (*
DAE_guess) (
double &,
double [],
double [],
double []);
double * return_float_space(void)
Gets the double array used by COLDAE.
int get_ODE_order(void)
Get the order of ODE's.
Definition: DAE_solver.h:209
double a
Leftmost boundary point.
Definition: DAE_solver.h:416
~DAESolver()
Destructor - clear all data.
Definition: DAE_solver.h:134
void(* gsub_ptr)(int &, double[], double &)
Definition: DAE_solver.h:48
void set_guess(void(*)(double &, double[], double[], double[]))
Sets the initial-guess function.
int zeta_size
Size of zeta array.
Definition: DAE_solver.h:341
void set_dfsub(void(*dfsub)(double &, double[], double[], double[]))
Sets the Jacobian of the ODE function.
double * tol
Array of tolerances.
Definition: DAE_solver.h:354
void set_tolerance(int ltol_size=0, int *ltol=NULL, double *tol=NULL)
Set tolerances for solution components.
int output_level
Output level.
Definition: DAE_solver.h:391
void(* dgsub_ptr)(int &, double[], double[])
Definition: DAE_solver.h:49
bool set_intialmeshsize
Set initial mesh flag.
Definition: DAE_solver.h:372
bool set_tol
tol flag
Definition: DAE_solver.h:356
void set_side_conditions(int zeta_size, double *zeta)
Sets additional boundary points.
bool have_dfsub
Boolian for jacobian of fsub wrapper function.
Definition: DAE_solver.h:329
void set_integer_space(int ispace_size=0)
Sets the size of the integer array used by COLDAE.
bool set_fspace
Set fspace flag.
Definition: DAE_solver.h:386
void(* DAE_guess)(double &, double[], double[], double[])
Pointer to geuss function.
Definition: DAE_solver.h:297
void set_collocation_points(int pnts)
Sets number of collocation points to be used in each subinterval.
int * ltol
Array used to hold location of tolerances.
Definition: DAE_solver.h:347
void for_to_c_matrix(int rows, int cols, double *fmat, double **cmat)
Converts a FORTRAN 2D array to a C/C++ 2D array.
void set_gsub(void(*)(int &, double[], double &))
Sets the boundary condition function.
int DAE_index
DAE index.
Definition: DAE_solver.h:406
bool set_ltol
Set ltol flag.
Definition: DAE_solver.h:349
void set_initial_mesh_size(int pnts)
Sets the initial mesh size.
int * return_integer_space(void)
Gets the integer array used by COLDAE.
void(* fsub_ptr)(double &, double[], double[], double[])
Definition: DAE_solver.h:46
void set_ipar(void)
Creates an array that contains various information about how to solve the DAE.
bool set_collpnts
Set collocation flag.
Definition: DAE_solver.h:367
int fspace_size
Size of fspace.
Definition: DAE_solver.h:384
int ispace_size
size of ispace
Definition: DAE_solver.h:377
void(* DAE_dgsub)(int &, double[], double[])
Definition: DAE_solver.h:295
bool have_guess
Boolian for guess function.
Definition: DAE_solver.h:335
void(* guess_ptr)(double &, double[], double[], double[])
Definition: DAE_solver.h:50
int fixpnt_size
Size of fixpnt array.
Definition: DAE_solver.h:398
bool have_fsub
Boolian for DAE function.
Definition: DAE_solver.h:327
DAESolver(int m_comp, int ny, int m[], double a, double b, double zeta[], void(*fsub)(double &, double[], double[], double[]), void(*dfsub)(double &, double[], double[], double[]), void(*gsub)(int &, double[], double &), void(*dgsub)(int &, double[], double[]), void(*guess)(double &, double[], double[], double[])=NULL)
Constructor.
void(* DAE_fsub)(double &, double[], double[], double[])
Definition: DAE_solver.h:289
double * zeta
Array of boundary locations.
Definition: DAE_solver.h:339
int ltol_size
Size of ltol array.
Definition: DAE_solver.h:351
int num_ODEs
Number of ODEs.
Definition: DAE_solver.h:409
This class provides an interface to the Fortran 77 code COLDAE.
Definition: DAE_solver.h:105
void coldae_(int &, int &, int[], double &, double &, double[], int[], int[], double[], double[], int[], double[], int &, void(*)(double &, double[], double[], double[]), void(*)(double &, double[], double[], double[]), void(*)(int &, double[], double &), void(*)(int &, double[], double[]), void(*)(double &, double[], double[], double[]))
bool use_cont
use simple continuation
Definition: DAE_solver.h:420
void(* dfsub_ptr)(double &, double[], double[], double[])
Definition: DAE_solver.h:47
void set_float_space(int fspace_size=0, double *fspace=NULL)
Sets the size and location of double array used by COLDAE.
void appsln_(double &, double[], double[], double[], int[])
void set_fixpnts(int fixpnt_size=1, double *fixpnt=NULL)
Set the fixed points in the mesh.
void set_linear(void)
Indicate if the problem is linear.
int solvercontrol
solver control parameter (see ipar(10))
Definition: DAE_solver.h:403
void DAE_dummy_guess(double &x, double z[], double y[], double df[])
A dummy guess function to be provided to COLDAE when a user wishes to provide none.
int collpnts
Number of collocation points.
Definition: DAE_solver.h:365
int * ODEs_Orders
A pointer to an array that contains the orders of each of the ODEs.
Definition: DAE_solver.h:414
int * ispace
Integer array used by COLDAE.
Definition: DAE_solver.h:375
int DAE_solve(void)
Solves the DAE.
void set_prob_size(int num_ODEs, int num_Alg_Const, int *ODEs_Orders)
Sets the problem size.
int * ipar
Array of info passed to COLDAE.
Definition: DAE_solver.h:359
void get_copy_final_mesh(double *mesh)
Returns a copy of the final mesh.
void use_sol_as_guess()
Overloads delete operator.
Definition: DAE_solver.h:258
int get_size_final_mesh(void)
Get the number of mesh points in the final mesh.
Definition: DAE_solver.h:245
void clear_mem(void)
Deletes dynamic memory.
void DAE_solution(double x, double z[], double y[])
Return a solution for the DAE at point x.
void c_to_for_matrix(int rows, int cols, double **cmat, double *fmat)
Converts a C/C++ 2D array to a Fortran 2D array.
double * fixpnt
Array of fixed points.
Definition: DAE_solver.h:394
bool use_old_fspace
USer old fspace for continuation.
Definition: DAE_solver.h:388
bool set_zeta
Set zeta flag.
Definition: DAE_solver.h:343
bool have_gsub
Boolian for boundary condition function.
Definition: DAE_solver.h:331
void set_dgsub(void(*)(int &, double[], double[]))
Sets the partial derivative of the boundary condition function.
void use_simple_cont(void)
Allows for simple continuation.
void set_output(int level)
Set output level for COLDAE.
void(* DAE_dfsub)(double &, double[], double[], double[])
Pointer to wrapper for jacobian of fsub.
Definition: DAE_solver.h:291
int intialmeshsize
Size of initial mesh.
Definition: DAE_solver.h:370
double * fspace
Double array used by COLDAE.
Definition: DAE_solver.h:382
bool have_dgsub
Boolian for jacobian of boundary condition function.
Definition: DAE_solver.h:333
bool set_ispace
ispace flag
Definition: DAE_solver.h:379
void set_fsub(void(*fsub)(double &, double[], double[], double[]))
Set the ODE function.
double b
Rightmost boundary point.
Definition: DAE_solver.h:418
void set_solver_control(int control)
Set the solver control parameter.
bool set_solvercontrol
solver control flag
Definition: DAE_solver.h:401
void set_boundary_points(int a, int b)
Sets the left and right-most boundary points.
void set_DAE_index(int index)
Set index of DAE.
void(* DAE_gsub)(int &, double[], double &)
Pointer to boundary condition function.
Definition: DAE_solver.h:293
bool linear
Linear flag.
Definition: DAE_solver.h:362
bool set_fixpnt
fixpnt flag
Definition: DAE_solver.h:396
int num_Alg_Const
Number of algebraic constants.
Definition: DAE_solver.h:411