DPsim
ODESolver.cpp
1 /* Copyright 2017-2021 Institute for Automation of Complex Power Systems,
2  * EONERC, RWTH Aachen University
3  *
4  * This Source Code Form is subject to the terms of the Mozilla Public
5  * License, v. 2.0. If a copy of the MPL was not distributed with this
6  * file, You can obtain one at https://mozilla.org/MPL/2.0/.
7  *********************************************************************************/
8 
9 #include <dpsim-models/SimPowerComp.h>
10 #include <dpsim/ODESolver.h>
11 
12 using namespace DPsim;
13 
14 ODESolver::ODESolver(String name, const CPS::ODEInterface::Ptr &comp,
15  bool implicit_integration, Real timestep)
16  : Solver(name, CPS::Logger::Level::info), mComponent(comp),
17  mImplicitIntegration(implicit_integration), mTimestep(timestep) {
18  mProbDim = mComponent->mOdePreState->get().rows();
19  initialize();
20 }
21 
23  mStates = N_VNew_Serial(mProbDim);
24  // Set initial value: (Different from DAESolver), only for already initialized components!
25  // XXX
26  N_VSetArrayPointer((**mComponent->mOdePostState).data(), mStates);
27  // Forbid SUNdials from deleting the underlying state vector (which is managed
28  // by our attribute / shared_ptr system)
29  NV_OWN_DATA_S(mStates) = false;
30 
31  // Analogous to DAESolver
32  CPS::ODEInterface::Ptr dummy = mComponent;
33  mStSpFunction = [dummy](double t, const double y[], double ydot[]) {
34  dummy->odeStateSpace(t, y, ydot);
35  };
36  mJacFunction = [dummy](double t, const double y[], double fy[], double J[],
37  double tmp1[], double tmp2[], double tmp3[]) {
38  dummy->odeJacobian(t, y, fy, J, tmp1, tmp2, tmp3);
39  };
40 
41  // Causes numerical issues, better allocate in every step-> see step
42  /*mArkode_mem= ARKodeCreate();
43  if (check_flag(mArkode_mem, "ARKodeCreate", 0))
44  mFlag=1;
45 
46  mFlag = ARKodeSetUserData(mArkode_mem, this);
47  if (check_flag(&mFlag, "ARKodeSetUserData", 1))
48  mFlag=1;*/
49 
50  /* Call ARKodeInit to initialize the integrator memory and specify the
51  right-hand side function in y'=f(t,y), the inital time T0, and
52  the initial dependent variable vector y(fluxes+mech. vars).*/
53  /* if(mImplicitIntegration){
54  mFlag = ARKodeInit(mArkode_mem, NULL, &ODESolver::StateSpaceWrapper, initial_time, mStates);
55  if (check_flag(&mFlag, "ARKodeInit", 1)) throw CPS::Exception();
56 
57  // Initialize dense matrix data structure
58  A = SUNDenseMatrix(mProbDim, mProbDim);
59  if (check_flag((void *)A, "SUNDenseMatrix", 0)) throw CPS::Exception();
60 
61  // Initialize linear solver
62  LS = SUNDenseLinearSolver(mStates, A);
63  if (check_flag((void *)LS, "SUNDenseLinearSolver", 0)) throw CPS::Exception();
64 
65  // Attach matrix and linear solver
66  mFlag = ARKDlsSetLinearSolver(mArkode_mem, LS, A);
67  if (check_flag(&mFlag, "ARKDlsSetLinearSolver", 1)) throw CPS::Exception();
68 
69  // Set Jacobian routine
70  mFlag = ARKDlsSetJacFn(mArkode_mem, &ODESolver::JacobianWrapper);
71  if (check_flag(&mFlag, "ARKDlsSetJacFn", 1)) throw CPS::Exception();
72  }
73  else {
74  mFlag = ARKodeInit(mArkode_mem, &ODESolver::StateSpaceWrapper, NULL, initial_time, mStates);
75  if (check_flag(&mFlag, "ARKodeInit", 1)) throw CPS::Exception();
76  }*/
77 
78  // Shifted to every step because of numerical issues
79 
80  // Specify Runge-Kutta Method/order
81  /*mFlag = ARKodeSetOrder(mArkode_mem, 4);
82  if (check_flag(&mFlag, "ARKodeOrderSet", 1))
83  mFlag=1;*/
84 
85  /*mFlag = ARKodeSStolerances(mArkode_mem, reltol, abstol);
86  if (check_flag(&mFlag, "ARKodeSStolerances", 1))
87  mFlag=1;*/
88 }
89 
90 int ODESolver::StateSpaceWrapper(realtype t, N_Vector y, N_Vector ydot,
91  void *user_data) {
92  ODESolver *self = reinterpret_cast<ODESolver *>(user_data);
93  return self->StateSpace(t, y, ydot);
94 }
95 
96 int ODESolver::StateSpace(realtype t, N_Vector y, N_Vector ydot) {
97  mStSpFunction(t, NV_DATA_S(y), NV_DATA_S(ydot));
98  return 0;
99 }
100 
101 int ODESolver::JacobianWrapper(realtype t, N_Vector y, N_Vector fy, SUNMatrix J,
102  void *user_data, N_Vector tmp1, N_Vector tmp2,
103  N_Vector tmp3) {
104  ODESolver *self = reinterpret_cast<ODESolver *>(user_data);
105  return self->Jacobian(t, y, fy, J, tmp1, tmp2, tmp3);
106 }
107 
108 int ODESolver::Jacobian(realtype t, N_Vector y, N_Vector fy, SUNMatrix J,
109  N_Vector tmp1, N_Vector tmp2, N_Vector tmp3) {
110  mJacFunction(t, NV_DATA_S(y), NV_DATA_S(fy), SM_DATA_D(J), NV_DATA_S(tmp1),
111  NV_DATA_S(tmp2), NV_DATA_S(tmp3));
112  return 0;
113 }
114 
115 Real ODESolver::step(Real initial_time) {
116  // Not absolutely necessary; realtype by default double (same as Real)
117  realtype T0 = (realtype)initial_time;
118  realtype Tf = (realtype)initial_time + mTimestep;
119 
121  long int nst;
123  long int netf;
124 
125  mComponent->mOdePostState->set(mComponent->mOdePreState->get());
126 
127  // Better allocate the arkode memory here to prevent numerical problems
128  mArkode_mem = ARKodeCreate();
129  if (check_flag(mArkode_mem, "ARKodeCreate", 0))
130  mFlag = 1;
131 
132  mFlag = ARKodeSetUserData(mArkode_mem, this);
133  if (check_flag(&mFlag, "ARKodeSetUserData", 1))
134  mFlag = 1;
135 
136  /* Call ARKodeInit to initialize the integrator memory and specify the
137  right-hand side function in y'=f(t,y), the inital time T0, and
138  the initial dependent variable vector y(fluxes+mech. vars).*/
139  if (mImplicitIntegration) {
140  mFlag = ARKodeInit(mArkode_mem, NULL, &ODESolver::StateSpaceWrapper,
141  initial_time, mStates);
142  if (check_flag(&mFlag, "ARKodeInit", 1))
143  throw CPS::Exception();
144 
145  // Initialize dense matrix data structure
146  A = SUNDenseMatrix(mProbDim, mProbDim);
147  if (check_flag((void *)A, "SUNDenseMatrix", 0))
148  throw CPS::Exception();
149 
150  // Initialize linear solver
151  LS = SUNDenseLinearSolver(mStates, A);
152  if (check_flag((void *)LS, "SUNDenseLinearSolver", 0))
153  throw CPS::Exception();
154 
155  // Attach matrix and linear solver
156  mFlag = ARKDlsSetLinearSolver(mArkode_mem, LS, A);
157  if (check_flag(&mFlag, "ARKDlsSetLinearSolver", 1))
158  throw CPS::Exception();
159 
160  // Set Jacobian routine
161  mFlag = ARKDlsSetJacFn(mArkode_mem, &ODESolver::JacobianWrapper);
162  if (check_flag(&mFlag, "ARKDlsSetJacFn", 1))
163  throw CPS::Exception();
164  } else {
165  mFlag = ARKodeInit(mArkode_mem, &ODESolver::StateSpaceWrapper, NULL,
166  initial_time, mStates);
167  if (check_flag(&mFlag, "ARKodeInit", 1))
168  throw CPS::Exception();
169  }
170 
171  mFlag = ARKodeSStolerances(mArkode_mem, reltol, abstol);
172  if (check_flag(&mFlag, "ARKodeSStolerances", 1))
173  mFlag = 1;
174 
175  // Main integrator loop
176  realtype t = T0;
177  while (Tf - t > 1.0e-15) {
178  mFlag = ARKode(mArkode_mem, Tf, mStates, &t, ARK_NORMAL);
179  if (check_flag(&mFlag, "ARKode", 1))
180  break;
181  }
182 
183  // Get some statistics to check for numerical problems (instability, blow-up etc)
184  mFlag = ARKodeGetNumSteps(mArkode_mem, &nst);
185  if (check_flag(&mFlag, "ARKodeGetNumSteps", 1))
186  return 1;
187  mFlag = ARKodeGetNumErrTestFails(mArkode_mem, &netf);
188  if (check_flag(&mFlag, "ARKodeGetNumErrTestFails", 1))
189  return 1;
190 
191  ARKodeFree(&mArkode_mem);
192  SUNLinSolFree(LS);
193  SUNMatDestroy(A);
194 
195  // Print statistics:
196  //std::cout << "Number Computing Steps: "<< nst << " Number Error-Test-Fails: " << netf << std::endl;
197  return Tf;
198 }
199 
200 void ODESolver::SolveTask::execute(Real time, Int timeStepCount) {
201  mSolver.step(time);
202 }
203 
204 // ARKode-Error checking functions
205 // Check function return value...
206 // opt == 0 means SUNDIALS function allocates memory so check if
207 // returned NULL pointer
208 // opt == 1 means SUNDIALS function returns a flag so check if
209 // flag >= 0
210 // opt == 2 means function allocates memory so check if returned
211 // NULL pointer
212 int ODESolver::check_flag(void *flagvalue, const std::string &funcname,
213  int opt) {
214  int *errflag;
215 
216  // Check if SUNDIALS function returned NULL pointer - no memory allocated
217  if (opt == 0 && flagvalue == NULL) {
218  std::cout << "\nSUNDIALS_ERROR: " << funcname
219  << " failed - returned NULL pointer\n\n";
220  return 1;
221  }
222 
223  // Check if flag < 0
224  else if (opt == 1) {
225  errflag = (int *)flagvalue;
226  if (*errflag < 0) {
227  std::cout << "\nSUNDIALS_ERROR: " << funcname
228  << " failed with flag = " << *errflag << "\n\n";
229  return 1;
230  }
231  }
232 
233  // Check if function returned NULL pointer - no memory allocated
234  else if (opt == 2 && flagvalue == NULL) {
235  std::cout << "\nMEMORY_ERROR: " << funcname
236  << " failed - returned NULL pointer\n\n";
237  return 1;
238  }
239 
240  return 0;
241 }
242 
243 ODESolver::~ODESolver() { N_VDestroy(mStates); }
Solver class for ODE (Ordinary Differential Equation) systems.
Definition: ODESolver.h:25
ODESolver(String name, const CPS::ODEInterface::Ptr &comp, bool implicit_integration, Real timestep)
Create solve object with corresponding component and information on the integration type.
Definition: ODESolver.cpp:14
CPS::ODEInterface::Ptr mComponent
Component to simulate, possible specialized component needed.
Definition: ODESolver.h:28
Real step(Real initial_time)
Solve system for the current time.
Definition: ODESolver.cpp:115
void * mArkode_mem
Memory block allocated by ARKode.
Definition: ODESolver.h:35
int mFlag
Template Jacobian Matrix (implicit solver)
Definition: ODESolver.h:63
realtype abstol
Scalar absolute tolerance.
Definition: ODESolver.h:54
SUNLinearSolver LS
Empty linear solver object.
Definition: ODESolver.h:45
int check_flag(void *flagvalue, const std::string &funcname, int opt)
ARKode- standard error detection function; in DAE-solver not detection function is used -> for effici...
Definition: ODESolver.cpp:212
N_Vector mStates
State vector.
Definition: ODESolver.h:37
void initialize()
Initialize ARKode-solve_environment.
Definition: ODESolver.cpp:22
bool mImplicitIntegration
Indicates whether the ODE shall be solved using an implicit scheme.
Definition: ODESolver.h:41
~ODESolver()
Deallocate all memory.
Definition: ODESolver.cpp:243
Real mTimestep
Constant time step.
Definition: ODESolver.h:48
SUNMatrix A
Empty matrix for linear solve in each Newton step while solving the Jacobian Matrix.
Definition: ODESolver.h:43
Int mProbDim
Number of differential Variables (states)
Definition: ODESolver.h:31
realtype reltol
Relative tolerance.
Definition: ODESolver.h:52
static int StateSpaceWrapper(realtype t, N_Vector y, N_Vector ydot, void *user_data)
use wrappers similar to DAE_Solver
Definition: ODESolver.cpp:90
Base class for more specific solvers such as MNA, ODE or IDA.
Definition: Solver.h:30