HyPar  1.0
Finite-Difference Hyperbolic-Parabolic PDE Solver on Cartesian Grids
FPPowerSystem3BusInitialize.c File Reference

Function to initialize the 3-bus power system model. More...

#include <stdio.h>
#include <stdlib.h>
#include <math.h>
#include <string.h>
#include <basic.h>
#include <arrayfunctions.h>
#include <physicalmodels/fppowersystem3bus.h>
#include <mpivars.h>
#include <hypar.h>

Go to the source code of this file.

Functions

double FPPowerSystem3BusComputeCFL (void *, void *, double, double)
 
double FPPowerSystem3BusComputeDiffNumber (void *, void *, double, double)
 
int FPPowerSystem3BusAdvection (double *, double *, int, void *, double)
 
int FPPowerSystem3BusDiffusion (double *, double *, int, int, void *, double)
 
int FPPowerSystem3BusUpwind (double *, double *, double *, double *, double *, double *, int, void *, double)
 
int FPPowerSystem3BusInitialize (void *s, void *m)
 

Detailed Description

Function to initialize the 3-bus power system model.

Author
Debojyoti Ghosh

Definition in file FPPowerSystem3BusInitialize.c.

Function Documentation

◆ FPPowerSystem3BusComputeCFL()

double FPPowerSystem3BusComputeCFL ( void *  s,
void *  m,
double  dt,
double  t 
)

Computes the maximum CFL number over the domain. Note that the CFL is computed over the local domain on this processor only.

Parameters
sSolver object of type HyPar
mMPI object of type MPIVariables
dtTime step size for which to compute the CFL
tTime

Definition at line 19 of file FPPowerSystem3BusComputeCFL.c.

25 {
26  HyPar *solver = (HyPar*) s;
27  FPPowerSystem3Bus *params = (FPPowerSystem3Bus*) solver->physics;
28 
29  int ndims = solver->ndims;
30  int ghosts = solver->ghosts;
31  int *dim = solver->dim_local;
32 
33  double max_cfl = 0;
34  int index[ndims];
35  int done = 0; _ArraySetValue_(index,ndims,0);
36  while (!done) {
37  double x[ndims], dxinv[ndims], drift[ndims];
38  _GetCoordinate_(0,index[0],dim,ghosts,solver->x,x[0]);
39  _GetCoordinate_(1,index[1],dim,ghosts,solver->x,x[1]);
40  _GetCoordinate_(2,index[2],dim,ghosts,solver->x,x[2]);
41  _GetCoordinate_(3,index[3],dim,ghosts,solver->x,x[3]);
42  _GetCoordinate_(0,index[0],dim,ghosts,solver->dxinv,dxinv[0]);
43  _GetCoordinate_(1,index[1],dim,ghosts,solver->dxinv,dxinv[1]);
44  _GetCoordinate_(2,index[2],dim,ghosts,solver->dxinv,dxinv[2]);
45  _GetCoordinate_(3,index[3],dim,ghosts,solver->dxinv,dxinv[3]);
46 
48 
49  double local_cfl[ndims];
50  local_cfl[0] = absolute(drift[0]) * dt * dxinv[0];
51  local_cfl[1] = absolute(drift[1]) * dt * dxinv[1];
52  local_cfl[2] = absolute(drift[2]) * dt * dxinv[2];
53  local_cfl[3] = absolute(drift[3]) * dt * dxinv[3];
54 
55  if (local_cfl[0] > max_cfl) max_cfl = local_cfl[0];
56  if (local_cfl[1] > max_cfl) max_cfl = local_cfl[1];
57  if (local_cfl[2] > max_cfl) max_cfl = local_cfl[2];
58  if (local_cfl[3] > max_cfl) max_cfl = local_cfl[3];
59 
60  _ArrayIncrementIndex_(ndims,dim,index,done);
61  }
62 
63  return(max_cfl);
64 }
#define absolute(a)
Definition: math_ops.h:32
double * x
Definition: hypar.h:107
#define drift(x)
Definition: fpdoublewell.h:31
int ndims
Definition: hypar.h:26
Structure containing all solver-specific variables and functions.
Definition: hypar.h:23
int FPPowerSystem3BusDriftFunction(int, void *, double *, double, double *)
#define _ArraySetValue_(x, size, value)
int * dim_local
Definition: hypar.h:37
#define _ArrayIncrementIndex_(N, imax, i, done)
#define _GetCoordinate_(dir, i, dim, ghosts, x, coord)
Definition: basic.h:31
void * physics
Definition: hypar.h:266
int ghosts
Definition: hypar.h:52
Structure containing variable and parameters specific to the 3-bus power system model. This structure contains the physical parameters and variables for the Fokker-Planck model for a 3-bus power system.
double * dxinv
Definition: hypar.h:110

◆ FPPowerSystem3BusComputeDiffNumber()

double FPPowerSystem3BusComputeDiffNumber ( void *  s,
void *  m,
double  dt,
double  t 
)

Computes the maximum diffusion number over the domain. Note that the diffusion is computed over the local domain on this processor only.

Parameters
sSolver object of type HyPar
mMPI object of type MPIVariables
dtTime step size for which to compute the CFL
tTime

Definition at line 19 of file FPPowerSystem3BusComputeDiffNumber.c.

25 {
26  HyPar *solver = (HyPar*) s;
27  FPPowerSystem3Bus *params = (FPPowerSystem3Bus*) solver->physics;
28 
29  int ndims = solver->ndims;
30  int ghosts = solver->ghosts;
31  int *dim = solver->dim_local;
32 
33  double max_diff = 0;
34  int index[ndims];
35  int done = 0; _ArraySetValue_(index,ndims,0);
36  while (!done) {
37  double dxinv[ndims],dissp[ndims*ndims];
38  _GetCoordinate_(0,index[0],dim,ghosts,solver->dxinv,dxinv[0]);
39  _GetCoordinate_(1,index[1],dim,ghosts,solver->dxinv,dxinv[1]);
40  _GetCoordinate_(2,index[2],dim,ghosts,solver->dxinv,dxinv[2]);
41  _GetCoordinate_(3,index[3],dim,ghosts,solver->dxinv,dxinv[3]);
42  FPPowerSystem3BusDissipationFunction(0,0,params,t,dissp);
43 
44  double local_diff[ndims];
45  local_diff[0] = absolute(dissp[0*ndims+0]) * dt * dxinv[0] * dxinv[0];
46  local_diff[1] = absolute(dissp[1*ndims+1]) * dt * dxinv[1] * dxinv[1];
47  local_diff[2] = absolute(dissp[2*ndims+2]) * dt * dxinv[2] * dxinv[2];
48  local_diff[3] = absolute(dissp[3*ndims+3]) * dt * dxinv[3] * dxinv[3];
49 
50  if (local_diff[0] > max_diff) max_diff = local_diff[0];
51  if (local_diff[1] > max_diff) max_diff = local_diff[1];
52  if (local_diff[2] > max_diff) max_diff = local_diff[2];
53  if (local_diff[3] > max_diff) max_diff = local_diff[3];
54 
55  _ArrayIncrementIndex_(ndims,dim,index,done);
56  }
57 
58  return(max_diff);
59 }
#define absolute(a)
Definition: math_ops.h:32
int FPPowerSystem3BusDissipationFunction(int, int, void *, double, double *)
int ndims
Definition: hypar.h:26
Structure containing all solver-specific variables and functions.
Definition: hypar.h:23
#define _ArraySetValue_(x, size, value)
int * dim_local
Definition: hypar.h:37
#define _ArrayIncrementIndex_(N, imax, i, done)
#define _GetCoordinate_(dir, i, dim, ghosts, x, coord)
Definition: basic.h:31
void * physics
Definition: hypar.h:266
int ghosts
Definition: hypar.h:52
Structure containing variable and parameters specific to the 3-bus power system model. This structure contains the physical parameters and variables for the Fokker-Planck model for a 3-bus power system.
double * dxinv
Definition: hypar.h:110

◆ FPPowerSystem3BusAdvection()

int FPPowerSystem3BusAdvection ( double *  f,
double *  u,
int  dir,
void *  s,
double  t 
)

Compute the advection term for the FPPowerSystem3Bus system: Since the advection coefficient is a function of x and not the solution, here the flux is set to the solution. The advection velocity is multiplied during upwinding FPPowerSystem3BusUpwind().

Parameters
fArray to hold the computed flux vector (same layout as u)
uArray with the solution vector
dirSpatial dimension for which to compute the flux
sSolver object of type HyPar
tCurrent simulation time

Definition at line 16 of file FPPowerSystem3BusAdvection.c.

23 {
24  HyPar *solver = (HyPar*) s;
26  return(0);
27 }
Structure containing all solver-specific variables and functions.
Definition: hypar.h:23
int npoints_local_wghosts
Definition: hypar.h:42
#define _ArrayCopy1D_(x, y, size)

◆ FPPowerSystem3BusDiffusion()

int FPPowerSystem3BusDiffusion ( double *  f,
double *  u,
int  dir1,
int  dir2,
void *  s,
double  t 
)

Compute the dissipation term for the FPPowerSystem3Bus system

Parameters
fArray to hold the computed dissipation term vector (same layout as u)
uArray with the solution vector
dir1First spatial dimension for the dissipation term being computed
dir2Second spatial dimension for the dissipation term being computed
sSolver object of type HyPar
tCurrent simulation time

Definition at line 15 of file FPPowerSystem3BusDiffusion.c.

23 {
24  HyPar *solver = (HyPar*) s;
25  FPPowerSystem3Bus *params = (FPPowerSystem3Bus*) solver->physics;
26  int i, v;
27 
28  int *dim = solver->dim_local;
29  int ghosts = solver->ghosts;
30  static int index[_MODEL_NDIMS_], bounds[_MODEL_NDIMS_], offset[_MODEL_NDIMS_];
31  static double dissipation[_MODEL_NDIMS_*_MODEL_NDIMS_];
32 
33  /* set bounds for array index to include ghost points */
34  _ArrayCopy1D_(dim,bounds,_MODEL_NDIMS_);
35  for (i=0; i<_MODEL_NDIMS_; i++) bounds[i] += 2*ghosts;
36 
37  /* set offset such that index is compatible with ghost point arrangement */
38  _ArraySetValue_(offset,_MODEL_NDIMS_,-ghosts);
39 
40  int done = 0; _ArraySetValue_(index,_MODEL_NDIMS_,0);
41  while (!done) {
42  int p; _ArrayIndex1DWO_(_MODEL_NDIMS_,dim,index,offset,ghosts,p);
43  FPPowerSystem3BusDissipationFunction(dir1,dir2,params,t,dissipation);
44  for (v = 0; v < _MODEL_NVARS_; v++) f[_MODEL_NVARS_*p+v] = dissipation[dir1*_MODEL_NDIMS_+dir2] * u[_MODEL_NVARS_*p+v];
45  _ArrayIncrementIndex_(_MODEL_NDIMS_,bounds,index,done);
46  }
47 
48  return(0);
49 }
Structure containing all solver-specific variables and functions.
Definition: hypar.h:23
#define _MODEL_NDIMS_
Definition: euler1d.h:56
int FPPowerSystem3BusDissipationFunction(int, int, void *, double, double *)
#define _ArrayIndex1DWO_(N, imax, i, offset, ghost, index)
#define _ArraySetValue_(x, size, value)
int * dim_local
Definition: hypar.h:37
#define _ArrayIncrementIndex_(N, imax, i, done)
void * physics
Definition: hypar.h:266
int ghosts
Definition: hypar.h:52
Structure containing variable and parameters specific to the 3-bus power system model. This structure contains the physical parameters and variables for the Fokker-Planck model for a 3-bus power system.
#define _MODEL_NVARS_
Definition: euler1d.h:58
#define _ArrayCopy1D_(x, y, size)

◆ FPPowerSystem3BusUpwind()

int FPPowerSystem3BusUpwind ( double *  fI,
double *  fL,
double *  fR,
double *  uL,
double *  uR,
double *  u,
int  dir,
void *  s,
double  t 
)

Compute the upwind flux at the interface, based on the drift velocity at that interface, from the left and right biased approximations to the interface flux. The drift (advection) velocity is multiplied to the solution in this function to get the advective flux.

Parameters
fIComputed upwind interface flux
fLLeft-biased reconstructed interface flux
fRRight-biased reconstructed interface flux
uLLeft-biased reconstructed interface solution
uRRight-biased reconstructed interface solution
uCell-centered solution
dirSpatial dimension (x or y)
sSolver object of type HyPar
tCurrent solution time

Definition at line 18 of file FPPowerSystem3BusUpwind.c.

29 {
30  HyPar *solver = (HyPar*) s;
31  FPPowerSystem3Bus *params = (FPPowerSystem3Bus*) solver->physics;
32  int done,v;
33 
34  int ndims = solver->ndims;
35  int nvars = solver->nvars;
36  int ghosts = solver->ghosts;
37  int *dim = solver->dim_local;
38 
39  int index_outer[ndims], index_inter[ndims], bounds_outer[ndims], bounds_inter[ndims];
40  _ArrayCopy1D_(dim,bounds_outer,ndims); bounds_outer[dir] = 1;
41  _ArrayCopy1D_(dim,bounds_inter,ndims); bounds_inter[dir] += 1;
42 
43  done = 0; _ArraySetValue_(index_outer,ndims,0);
44  while (!done) {
45  _ArrayCopy1D_(index_outer,index_inter,ndims);
46  for (index_inter[dir] = 0; index_inter[dir] < bounds_inter[dir]; index_inter[dir]++) {
47  int p; _ArrayIndex1D_(ndims,bounds_inter,index_inter,0, p);
48  double x[ndims]; /* coordinates of the interface */
49  double x1, x2, drift[ndims];
50  if (dir == 0) {
51  _GetCoordinate_(0,index_inter[0]-1,dim,ghosts,solver->x,x1);
52  _GetCoordinate_(0,index_inter[0] ,dim,ghosts,solver->x,x2);
53  x[0] = 0.5 * ( x1 + x2 );
54  _GetCoordinate_(1,index_inter[1],dim,ghosts,solver->x,x[1]);
55  _GetCoordinate_(2,index_inter[2],dim,ghosts,solver->x,x[2]);
56  _GetCoordinate_(3,index_inter[3],dim,ghosts,solver->x,x[3]);
57  } else if (dir == 1) {
58  _GetCoordinate_(0,index_inter[0],dim,ghosts,solver->x,x[0]);
59  _GetCoordinate_(1,index_inter[1]-1,dim,ghosts,solver->x,x1);
60  _GetCoordinate_(1,index_inter[1] ,dim,ghosts,solver->x,x2);
61  x[1] = 0.5 * ( x1 + x2 );
62  _GetCoordinate_(2,index_inter[2],dim,ghosts,solver->x,x[2]);
63  _GetCoordinate_(3,index_inter[3],dim,ghosts,solver->x,x[3]);
64  } else if (dir == 2) {
65  _GetCoordinate_(0,index_inter[0],dim,ghosts,solver->x,x[0]);
66  _GetCoordinate_(1,index_inter[1],dim,ghosts,solver->x,x[1]);
67  _GetCoordinate_(2,index_inter[2]-1,dim,ghosts,solver->x,x1);
68  _GetCoordinate_(2,index_inter[2] ,dim,ghosts,solver->x,x2);
69  x[2] = 0.5 * ( x1 + x2 );
70  _GetCoordinate_(3,index_inter[3],dim,ghosts,solver->x,x[3]);
71  } else if (dir == 3) {
72  _GetCoordinate_(0,index_inter[0],dim,ghosts,solver->x,x[0]);
73  _GetCoordinate_(1,index_inter[1],dim,ghosts,solver->x,x[1]);
74  _GetCoordinate_(2,index_inter[2],dim,ghosts,solver->x,x[2]);
75  _GetCoordinate_(3,index_inter[3]-1,dim,ghosts,solver->x,x1);
76  _GetCoordinate_(3,index_inter[3] ,dim,ghosts,solver->x,x2);
77  x[3] = 0.5 * ( x1 + x2 );
78  }
79  FPPowerSystem3BusDriftFunction(dir,params,x,t,drift);
80  for (v = 0; v < nvars; v++)
81  fI[nvars*p+v] = drift[dir] * (drift[dir] > 0 ? fL[nvars*p+v] : fR[nvars*p+v] );
82  }
83  _ArrayIncrementIndex_(ndims,bounds_outer,index_outer,done);
84  }
85 
86  return(0);
87 }
int nvars
Definition: hypar.h:29
double * x
Definition: hypar.h:107
#define drift(x)
Definition: fpdoublewell.h:31
int ndims
Definition: hypar.h:26
Structure containing all solver-specific variables and functions.
Definition: hypar.h:23
int FPPowerSystem3BusDriftFunction(int, void *, double *, double, double *)
#define _ArrayIndex1D_(N, imax, i, ghost, index)
#define _ArraySetValue_(x, size, value)
int * dim_local
Definition: hypar.h:37
#define _ArrayIncrementIndex_(N, imax, i, done)
#define _GetCoordinate_(dir, i, dim, ghosts, x, coord)
Definition: basic.h:31
void * physics
Definition: hypar.h:266
int ghosts
Definition: hypar.h:52
Structure containing variable and parameters specific to the 3-bus power system model. This structure contains the physical parameters and variables for the Fokker-Planck model for a 3-bus power system.
#define _ArrayCopy1D_(x, y, size)

◆ FPPowerSystem3BusInitialize()

int FPPowerSystem3BusInitialize ( void *  s,
void *  m 
)

Initialize the 3-bus power system model: Sets the default parameters, read in and set physics-related parameters, and set the physics-related function pointers in HyPar.

Parameters
sSolver object of type HyPar
mMPI object of type MPIVariables

Definition at line 27 of file FPPowerSystem3BusInitialize.c.

31 {
32  HyPar *solver = (HyPar*) s;
33  MPIVariables *mpi = (MPIVariables*) m;
34  FPPowerSystem3Bus *physics = (FPPowerSystem3Bus*) solver->physics;
35  int ferr, N;
37 
38  static int count = 0;
39 
40  if (solver->nvars != _MODEL_NVARS_) {
41  fprintf(stderr,"Error in FPPowerSystem3BusInitialize(): nvars has to be %d.\n",_MODEL_NVARS_);
42  return(1);
43  }
44  if (solver->ndims != _MODEL_NDIMS_) {
45  fprintf(stderr,"Error in FPPowerSystem3BusInitialize(): ndims has to be %d.\n",_MODEL_NDIMS_);
46  return(1);
47  }
48 
49  double pi = 4.0*atan(1.0);
50  /* default values of model parameters */
51  physics->N = 2;
52  physics->Pm1_avg = 0.8;
53  physics->Pm2_avg = 1.6;
54  physics->Pmref_avg = 0.79330781761651;
55  physics->H1 = 3.20;
56  physics->H2 = 6.40;
57  physics->Href = 13.60;
58  physics->E1 = 1.01556070860155;
59  physics->E2 = 1.0491099265981;
60  physics->Eref = 1.05623172878954;
61  physics->omegaB = 2*pi*60.0;
62  physics->sigma[0][0] = 0.0125;
63  physics->sigma[0][1] = 0.0;
64  physics->sigma[1][0] = 0.0;
65  physics->sigma[1][1] = 0.0125;
66  physics->lambda[0][0] = 10.0/physics->omegaB;
67  physics->lambda[0][1] = 0.0;
68  physics->lambda[1][0] = 0.0;
69  physics->lambda[1][1] = 10.0/physics->omegaB;
70  physics->gamma = 0.25;
71 
72  physics->G = (double*) calloc (3*3,sizeof(double));
73  physics->B = (double*) calloc (3*3,sizeof(double));
74 
75  physics->G[0*3+0] = 0.276805493111691;
76  physics->G[0*3+1] = 0.213024867595501;
77  physics->G[0*3+2] = 0.209205876527443;
78  physics->G[1*3+0] = 0.213024867595501;
79  physics->G[1*3+1] = 0.419642083051144;
80  physics->G[1*3+2] = 0.286592141665043;
81  physics->G[2*3+0] = 0.209205876527443;
82  physics->G[2*3+1] = 0.286592141665044;
83  physics->G[2*3+2] = 0.844559256324453;
84 
85  physics->B[0*3+0] = -2.36794416971567;
86  physics->B[0*3+1] = 1.08817493992579;
87  physics->B[0*3+2] = 1.22601259339234;
88  physics->B[1*3+0] = 1.08817493992579;
89  physics->B[1*3+1] = -2.72352378723346;
90  physics->B[1*3+2] = 1.51348094527252;
91  physics->B[2*3+0] = 1.22601259339234;
92  physics->B[2*3+1] = 1.51348094527252;
93  physics->B[2*3+2] = -2.98729895217208;
94 
95  /* reading physical model specific inputs - all processes */
96  FILE *in;
97  if ((!mpi->rank) && (!count)) printf("Reading physical model inputs from file \"physics.inp\".\n");
98  in = fopen("physics.inp","r");
99  if (!in) {
100  if (!mpi->rank) fprintf(stderr,"Error: File \"physics.inp\" not found.\n");
101  return(1);
102  } else {
103  char word[_MAX_STRING_SIZE_];
104  ferr = fscanf(in,"%s",word); if (ferr != 1) return(1);
105  if (!strcmp(word, "begin")){
106  while (strcmp(word, "end")){
107  ferr = fscanf(in,"%s",word); if (ferr != 1) return(1);
108  if (!strcmp(word,"Pm1_avg" )) {ferr=fscanf(in,"%lf",&physics->Pm1_avg ) ;if(ferr!=1)return(1);}
109  else if (!strcmp(word,"Pm2_avg" )) {ferr=fscanf(in,"%lf",&physics->Pm2_avg ) ;if(ferr!=1)return(1);}
110  else if (!strcmp(word,"Pmref_avg" )) {ferr=fscanf(in,"%lf",&physics->Pmref_avg ) ;if(ferr!=1)return(1);}
111  else if (!strcmp(word,"H1" )) {ferr=fscanf(in,"%lf",&physics->H1 ) ;if(ferr!=1)return(1);}
112  else if (!strcmp(word,"H2" )) {ferr=fscanf(in,"%lf",&physics->H2 ) ;if(ferr!=1)return(1);}
113  else if (!strcmp(word,"Href" )) {ferr=fscanf(in,"%lf",&physics->Href ) ;if(ferr!=1)return(1);}
114  else if (!strcmp(word,"E1" )) {ferr=fscanf(in,"%lf",&physics->E1 ) ;if(ferr!=1)return(1);}
115  else if (!strcmp(word,"E2" )) {ferr=fscanf(in,"%lf",&physics->E2 ) ;if(ferr!=1)return(1);}
116  else if (!strcmp(word,"Eref" )) {ferr=fscanf(in,"%lf",&physics->Eref ) ;if(ferr!=1)return(1);}
117  else if (!strcmp(word,"omegaB" )) {ferr=fscanf(in,"%lf",&physics->omegaB ) ;if(ferr!=1)return(1);}
118  else if (!strcmp(word,"gamma" )) {ferr=fscanf(in,"%lf",&physics->gamma ) ;if(ferr!=1)return(1);}
119  else if (!strcmp(word,"sigma" )) {
120  ferr=fscanf(in,"%lf",&physics->sigma[0][0]) ;if(ferr!=1)return(1);
121  ferr=fscanf(in,"%lf",&physics->sigma[0][1]) ;if(ferr!=1)return(1);
122  ferr=fscanf(in,"%lf",&physics->sigma[1][0]) ;if(ferr!=1)return(1);
123  ferr=fscanf(in,"%lf",&physics->sigma[1][1]) ;if(ferr!=1)return(1);
124  } else if (!strcmp(word,"lambda")) {
125  ferr=fscanf(in,"%lf",&physics->lambda[0][0]) ;if(ferr!=1)return(1);
126  ferr=fscanf(in,"%lf",&physics->lambda[0][1]) ;if(ferr!=1)return(1);
127  ferr=fscanf(in,"%lf",&physics->lambda[1][0]) ;if(ferr!=1)return(1);
128  ferr=fscanf(in,"%lf",&physics->lambda[1][1]) ;if(ferr!=1)return(1);
129  } else if (!strcmp(word,"G")) {
130  ferr=fscanf(in,"%lf",&physics->G[0*3+0]) ;if(ferr!=1)return(1);
131  ferr=fscanf(in,"%lf",&physics->G[0*3+1]) ;if(ferr!=1)return(1);
132  ferr=fscanf(in,"%lf",&physics->G[0*3+2]) ;if(ferr!=1)return(1);
133  ferr=fscanf(in,"%lf",&physics->G[1*3+0]) ;if(ferr!=1)return(1);
134  ferr=fscanf(in,"%lf",&physics->G[1*3+1]) ;if(ferr!=1)return(1);
135  ferr=fscanf(in,"%lf",&physics->G[1*3+2]) ;if(ferr!=1)return(1);
136  ferr=fscanf(in,"%lf",&physics->G[2*3+0]) ;if(ferr!=1)return(1);
137  ferr=fscanf(in,"%lf",&physics->G[2*3+1]) ;if(ferr!=1)return(1);
138  ferr=fscanf(in,"%lf",&physics->G[2*3+2]) ;if(ferr!=1)return(1);
139  } else if (!strcmp(word,"B")) {
140  ferr=fscanf(in,"%lf",&physics->B[0*3+0]) ;if(ferr!=1)return(1);
141  ferr=fscanf(in,"%lf",&physics->B[0*3+1]) ;if(ferr!=1)return(1);
142  ferr=fscanf(in,"%lf",&physics->B[0*3+2]) ;if(ferr!=1)return(1);
143  ferr=fscanf(in,"%lf",&physics->B[1*3+0]) ;if(ferr!=1)return(1);
144  ferr=fscanf(in,"%lf",&physics->B[1*3+1]) ;if(ferr!=1)return(1);
145  ferr=fscanf(in,"%lf",&physics->B[1*3+2]) ;if(ferr!=1)return(1);
146  ferr=fscanf(in,"%lf",&physics->B[2*3+0]) ;if(ferr!=1)return(1);
147  ferr=fscanf(in,"%lf",&physics->B[2*3+1]) ;if(ferr!=1)return(1);
148  ferr=fscanf(in,"%lf",&physics->B[2*3+2]) ;if(ferr!=1)return(1);
149  }
150  }
151  } else {
152  if (!mpi->rank) fprintf(stderr,"Error: Illegal format in file \"physics.inp\".\n");
153  return(1);
154  }
155  }
156  fclose(in);
157 
158  if (!strcmp(solver->SplitHyperbolicFlux,"yes")) {
159  if (!mpi->rank) {
160  fprintf(stderr,"Error in FPPowerSystem3BusInitialize: This physical model does not have a splitting ");
161  fprintf(stderr,"of the hyperbolic term defined.\n");
162  }
163  return(1);
164  }
165 
166  /* initializing physical model-specific functions */
172 
173  count++;
174  return(0);
175 }
double FPPowerSystem3BusComputeDiffNumber(void *, void *, double, double)
int nvars
Definition: hypar.h:29
char SplitHyperbolicFlux[_MAX_STRING_SIZE_]
Definition: hypar.h:92
int(* HFunction)(double *, double *, int, int, void *, double)
Definition: hypar.h:313
int FPPowerSystem3BusDiffusion(double *, double *, int, int, void *, double)
int(* Upwind)(double *, double *, double *, double *, double *, double *, int, void *, double)
Definition: hypar.h:295
int ndims
Definition: hypar.h:26
int FPPowerSystem3BusAdvection(double *, double *, int, void *, double)
int(* FFunction)(double *, double *, int, void *, double)
Definition: hypar.h:276
double(* ComputeCFL)(void *, void *, double, double)
Definition: hypar.h:269
Structure containing all solver-specific variables and functions.
Definition: hypar.h:23
double FPPowerSystem3BusComputeCFL(void *, void *, double, double)
#define _MAX_STRING_SIZE_
Definition: basic.h:14
#define _MODEL_NDIMS_
Definition: euler1d.h:56
int FPPowerSystem3BusUpwind(double *, double *, double *, double *, double *, double *, int, void *, double)
double(* ComputeDiffNumber)(void *, void *, double, double)
Definition: hypar.h:272
void * physics
Definition: hypar.h:266
Structure of MPI-related variables.
Structure containing variable and parameters specific to the 3-bus power system model. This structure contains the physical parameters and variables for the Fokker-Planck model for a 3-bus power system.
#define _MODEL_NVARS_
Definition: euler1d.h:58
#define _DECLARE_IERR_
Definition: basic.h:17