HyPar  1.0
Finite-Difference Hyperbolic-Parabolic PDE Solver on Cartesian Grids
FPPowerSystemInitialize.c File Reference
#include <stdio.h>
#include <stdlib.h>
#include <math.h>
#include <string.h>
#include <basic.h>
#include <arrayfunctions.h>
#include <physicalmodels/fppowersystem.h>
#include <mpivars.h>
#include <hypar.h>

Go to the source code of this file.

Functions

double FPPowerSystemComputeCFL (void *, void *, double, double)
 
double FPPowerSystemComputeDiffNumber (void *, void *, double, double)
 
int FPPowerSystemAdvection (double *, double *, int, void *, double)
 
int FPPowerSystemDiffusion (double *, double *, int, void *, double)
 
int FPPowerSystemUpwind (double *, double *, double *, double *, double *, double *, int, void *, double)
 
int FPPowerSystemPostStep (double *, void *, void *, double, int)
 
int FPPowerSystemPrintStep (void *, void *, double)
 
int FPPowerSystemInitialize (void *s, void *m)
 

Function Documentation

◆ FPPowerSystemComputeCFL()

double FPPowerSystemComputeCFL ( void *  ,
void *  ,
double  ,
double   
)

Definition at line 11 of file FPPowerSystemComputeCFL.c.

12 {
13  HyPar *solver = (HyPar*) s;
14  FPPowerSystem *params = (FPPowerSystem*)solver->physics;
15 
16  int ndims = solver->ndims;
17  int ghosts = solver->ghosts;
18  int *dim = solver->dim_local;
19 
20  double max_cfl = 0;
21  int index[ndims];
22  int done = 0; _ArraySetValue_(index,ndims,0);
23  while (!done) {
24  double x; _GetCoordinate_(0,index[0],dim,ghosts,solver->x,x);
25  double y; _GetCoordinate_(1,index[1],dim,ghosts,solver->x,y);
26  double dxinv; _GetCoordinate_(0,index[0],dim,ghosts,solver->dxinv,dxinv);
27  double dyinv; _GetCoordinate_(1,index[1],dim,ghosts,solver->dxinv,dyinv);
28  double drift_x= FPPowerSystemDriftFunction(0,params,x,y,t);
29  double drift_y= FPPowerSystemDriftFunction(1,params,x,y,t);
30 
31  double local_cfl_x = absolute(drift_x) * dt * dxinv;
32  double local_cfl_y = absolute(drift_y) * dt * dyinv;
33 
34  if (local_cfl_x > max_cfl) max_cfl = local_cfl_x;
35  if (local_cfl_y > max_cfl) max_cfl = local_cfl_y;
36 
37  _ArrayIncrementIndex_(ndims,dim,index,done);
38  }
39 
40  return(max_cfl);
41 }
#define absolute(a)
Definition: math_ops.h:32
double FPPowerSystemDriftFunction(int, void *, double, double, double)
double * x
Definition: hypar.h:107
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
double * dxinv
Definition: hypar.h:110

◆ FPPowerSystemComputeDiffNumber()

double FPPowerSystemComputeDiffNumber ( void *  ,
void *  ,
double  ,
double   
)

Definition at line 11 of file FPPowerSystemComputeDiffNumber.c.

12 {
13  HyPar *solver = (HyPar*) s;
14  FPPowerSystem *params = (FPPowerSystem*)solver->physics;
15 
16  int ndims = solver->ndims;
17  int ghosts = solver->ghosts;
18  int *dim = solver->dim_local;
19 
20  double max_diff = 0;
21  int index[ndims];
22  int done = 0; _ArraySetValue_(index,ndims,0);
23  while (!done) {
24  double dxinv; _GetCoordinate_(0,index[0],dim,ghosts,solver->dxinv,dxinv);
25  double dyinv; _GetCoordinate_(1,index[1],dim,ghosts,solver->dxinv,dyinv);
26  double dissp_x= FPPowerSystemDissipationFunction(0,params,t);
27  double dissp_y= FPPowerSystemDissipationFunction(1,params,t);
28 
29  double local_diff_x = absolute(dissp_x) * dt * dxinv * dxinv;
30  double local_diff_y = absolute(dissp_y) * dt * dyinv * dyinv;
31 
32  if (local_diff_x > max_diff) max_diff = local_diff_x;
33  if (local_diff_y > max_diff) max_diff = local_diff_y;
34 
35  _ArrayIncrementIndex_(ndims,dim,index,done);
36  }
37 
38  return(max_diff);
39 }
#define absolute(a)
Definition: math_ops.h:32
int ndims
Definition: hypar.h:26
double FPPowerSystemDissipationFunction(int, void *, double)
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
double * dxinv
Definition: hypar.h:110

◆ FPPowerSystemAdvection()

int FPPowerSystemAdvection ( double *  ,
double *  ,
int  ,
void *  ,
double   
)

Definition at line 7 of file FPPowerSystemAdvection.c.

8 {
9  HyPar *solver = (HyPar*) s;
10 // FPPowerSystem *params = (FPPowerSystem*)solver->physics;
11  int i, v;
12 
13  int *dim = solver->dim_local;
14  int ghosts = solver->ghosts;
15  int ndims = solver->ndims;
16  int nvars = solver->nvars;
17  int index[ndims], bounds[ndims], offset[ndims];
18 
19  /* set bounds for array index to include ghost points */
20  _ArrayCopy1D_(dim,bounds,ndims);
21  for (i=0; i<ndims; i++) bounds[i] += 2*ghosts;
22 
23  /* set offset such that index is compatible with ghost point arrangement */
24  _ArraySetValue_(offset,ndims,-ghosts);
25 
26  int done = 0; _ArraySetValue_(index,ndims,0);
27  while (!done) {
28  int p; _ArrayIndex1DWO_(ndims,dim,index,offset,ghosts,p);
29  for (v = 0; v < nvars; v++) f[nvars*p+v] = u[nvars*p+v];
30  _ArrayIncrementIndex_(ndims,bounds,index,done);
31  }
32 
33  return(0);
34 }
int nvars
Definition: hypar.h:29
int ndims
Definition: hypar.h:26
Structure containing all solver-specific variables and functions.
Definition: hypar.h:23
#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)
int ghosts
Definition: hypar.h:52
#define _ArrayCopy1D_(x, y, size)

◆ FPPowerSystemDiffusion()

int FPPowerSystemDiffusion ( double *  ,
double *  ,
int  ,
void *  ,
double   
)

Definition at line 9 of file FPPowerSystemDiffusion.c.

10 {
11  HyPar *solver = (HyPar*) s;
12  FPPowerSystem *params = (FPPowerSystem*) solver->physics;
13  int i, v;
14 
15  int *dim = solver->dim_local;
16  int ghosts = solver->ghosts;
17  int ndims = solver->ndims;
18  int nvars = solver->nvars;
19  int index[ndims], bounds[ndims], offset[ndims];
20 
21  /* set bounds for array index to include ghost points */
22  _ArrayCopy1D_(dim,bounds,ndims);
23  for (i=0; i<ndims; i++) bounds[i] += 2*ghosts;
24 
25  /* set offset such that index is compatible with ghost point arrangement */
26  _ArraySetValue_(offset,ndims,-ghosts);
27 
28  int done = 0; _ArraySetValue_(index,ndims,0);
29  while (!done) {
30  int p; _ArrayIndex1DWO_(ndims,dim,index,offset,ghosts,p);
31  double dissipation = FPPowerSystemDissipationFunction(dir,params,t);
32  for (v = 0; v < nvars; v++) f[nvars*p+v] = dissipation * u[nvars*p+v];
33  _ArrayIncrementIndex_(ndims,bounds,index,done);
34  }
35 
36  return(0);
37 }
int nvars
Definition: hypar.h:29
int ndims
Definition: hypar.h:26
double FPPowerSystemDissipationFunction(int, void *, double)
Structure containing all solver-specific variables and functions.
Definition: hypar.h:23
#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
#define _ArrayCopy1D_(x, y, size)

◆ FPPowerSystemUpwind()

int FPPowerSystemUpwind ( double *  ,
double *  ,
double *  ,
double *  ,
double *  ,
double *  ,
int  ,
void *  ,
double   
)

Definition at line 9 of file FPPowerSystemUpwind.c.

11 {
12  HyPar *solver = (HyPar*) s;
13  FPPowerSystem *params = (FPPowerSystem*)solver->physics;
14  int done,v;
15 
16  int ndims = solver->ndims;
17  int nvars = solver->nvars;
18  int ghosts = solver->ghosts;
19  int *dim = solver->dim_local;
20 
21  int index_outer[ndims], index_inter[ndims], bounds_outer[ndims], bounds_inter[ndims];
22  _ArrayCopy1D_(dim,bounds_outer,ndims); bounds_outer[dir] = 1;
23  _ArrayCopy1D_(dim,bounds_inter,ndims); bounds_inter[dir] += 1;
24 
25  done = 0; _ArraySetValue_(index_outer,ndims,0);
26  while (!done) {
27  _ArrayCopy1D_(index_outer,index_inter,ndims);
28  for (index_inter[dir] = 0; index_inter[dir] < bounds_inter[dir]; index_inter[dir]++) {
29  int p; _ArrayIndex1D_(ndims,bounds_inter,index_inter,0, p);
30  double x = 0,y = 0; /* x,y coordinates of the interface */
31  if (dir == 0) {
32  /* x-interface */
33  double x1, x2;
34  _GetCoordinate_(0,index_inter[0]-1,dim,ghosts,solver->x,x1);
35  _GetCoordinate_(0,index_inter[0] ,dim,ghosts,solver->x,x2);
36  x = 0.5 * ( x1 + x2 );
37  _GetCoordinate_(1,index_inter[1],dim,ghosts,solver->x,y);
38  } else if (dir == 1) {
39  /* y-interface */
40  _GetCoordinate_(0,index_inter[0],dim,ghosts,solver->x,x);
41  double y1, y2;
42  _GetCoordinate_(1,index_inter[1]-1,dim,ghosts,solver->x,y1);
43  _GetCoordinate_(1,index_inter[1] ,dim,ghosts,solver->x,y2);
44  y = 0.5 * ( y1 + y2 );
45  }
46  double drift = FPPowerSystemDriftFunction(dir,params,x,y,t);
47  for (v = 0; v < nvars; v++)
48  fI[nvars*p+v] = drift * (drift > 0 ? fL[nvars*p+v] : fR[nvars*p+v] );
49  }
50  _ArrayIncrementIndex_(ndims,bounds_outer,index_outer,done);
51  }
52 
53  return(0);
54 }
double FPPowerSystemDriftFunction(int, void *, double, double, double)
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
#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
#define _ArrayCopy1D_(x, y, size)

◆ FPPowerSystemPostStep()

int FPPowerSystemPostStep ( double *  ,
void *  ,
void *  ,
double  ,
int   
)

Definition at line 8 of file FPPowerSystemPostStep.c.

9 {
10  HyPar *solver = (HyPar*) s;
11  MPIVariables *mpi = (MPIVariables*) m;
12  FPPowerSystem *params = (FPPowerSystem*) solver->physics;
14 
15  int *dim = solver->dim_local;
16  int ghosts = solver->ghosts;
17  int ndims = solver->ndims;
18  int index[ndims];
19 
20  double local_sum = 0;
21  int done = 0; _ArraySetValue_(index,ndims,0);
22  while (!done) {
23  int p; _ArrayIndex1D_(ndims,dim,index,ghosts,p);
24  double dxinv; _GetCoordinate_(0,index[0],dim,ghosts,solver->dxinv,dxinv);
25  double dyinv; _GetCoordinate_(1,index[1],dim,ghosts,solver->dxinv,dyinv);
26  local_sum += (u[p] / (dxinv * dyinv));
27  _ArrayIncrementIndex_(ndims,dim,index,done);
28  }
29  double local_integral = local_sum;
30  double global_integral = 0;
31  IERR MPISum_double(&global_integral,&local_integral,1,&mpi->world); CHECKERR(ierr);
32  params->pdf_integral = global_integral;
33 
34  return(0);
35 }
#define IERR
Definition: basic.h:16
#define CHECKERR(ierr)
Definition: basic.h:18
int MPISum_double(double *, double *, int, void *)
Definition: MPISum.c:39
int ndims
Definition: hypar.h:26
double pdf_integral
Definition: fppowersystem.h:69
Structure containing all solver-specific variables and functions.
Definition: hypar.h:23
MPI_Comm world
#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 of MPI-related variables.
#define _DECLARE_IERR_
Definition: basic.h:17
double * dxinv
Definition: hypar.h:110

◆ FPPowerSystemPrintStep()

int FPPowerSystemPrintStep ( void *  ,
void *  ,
double   
)

Definition at line 5 of file FPPowerSystemPrintStep.c.

6 {
7  HyPar *solver = (HyPar*) s;
8  FPPowerSystem *params = (FPPowerSystem*) solver->physics;
9  printf("Domain integral of the probability density function: %1.16E\n",params->pdf_integral);
10  return(0);
11 }
double pdf_integral
Definition: fppowersystem.h:69
Structure containing all solver-specific variables and functions.
Definition: hypar.h:23
void * physics
Definition: hypar.h:266

◆ FPPowerSystemInitialize()

int FPPowerSystemInitialize ( void *  s,
void *  m 
)

Definition at line 20 of file FPPowerSystemInitialize.c.

21 {
22  HyPar *solver = (HyPar*) s;
23  MPIVariables *mpi = (MPIVariables*) m;
24  FPPowerSystem *physics = (FPPowerSystem*) solver->physics;
25  int ferr;
27 
28  static int count = 0;
29 
30  if (solver->nvars != _MODEL_NVARS_) {
31  fprintf(stderr,"Error in FPPowerSystemInitialize(): nvars has to be %d.\n",_MODEL_NVARS_);
32  return(1);
33  }
34  if (solver->ndims != _MODEL_NDIMS_) {
35  fprintf(stderr,"Error in FPPowerSystemInitialize(): ndims has to be %d.\n",_MODEL_NDIMS_);
36  return(1);
37  }
38 
39  /* default values of model parameters */
40  physics->H = 5.0;
41  physics->D = 5.0;
42  physics->E = 1.1358;
43  physics->V = 1.0;
44  physics->g1 = 0.545;
45  physics->g2 = 0.745;
46  physics->Pm = 0.9;
47  physics->l = 0.1;
48  physics->q = 1.0;
49  physics->O_s = 120*(4.0*atan(1.0));
50  physics->tf = 0.1;
51  physics->tcl = 0.2;
52 
53  /* reading physical model specific inputs - all processes */
54  FILE *in;
55  if ((!mpi->rank) && (!count)) printf("Reading physical model inputs from file \"physics.inp\".\n");
56  in = fopen("physics.inp","r");
57  if (!in) {
58  fprintf(stderr,"Error: File \"physics.inp\" not found.\n");
59  return(1);
60  } else {
61  char word[_MAX_STRING_SIZE_];
62  ferr = fscanf(in,"%s",word); if (ferr != 1) return(1);
63  if (!strcmp(word, "begin")){
64  while (strcmp(word, "end")){
65  ferr = fscanf(in,"%s",word); if (ferr != 1) return(1);
66  if (!strcmp(word, "inertia")) {ferr=fscanf(in,"%lf",&physics->H );if(ferr!=1)return(1);}
67  else if (!strcmp(word, "omega_s")) {ferr=fscanf(in,"%lf",&physics->O_s);if(ferr!=1)return(1);}
68  else if (!strcmp(word, "E" )) {ferr=fscanf(in,"%lf",&physics->E );if(ferr!=1)return(1);}
69  else if (!strcmp(word, "V" )) {ferr=fscanf(in,"%lf",&physics->V );if(ferr!=1)return(1);}
70  else if (!strcmp(word, "g1" )) {ferr=fscanf(in,"%lf",&physics->g1 );if(ferr!=1)return(1);}
71  else if (!strcmp(word, "g2" )) {ferr=fscanf(in,"%lf",&physics->g2 );if(ferr!=1)return(1);}
72  else if (!strcmp(word, "D" )) {ferr=fscanf(in,"%lf",&physics->D );if(ferr!=1)return(1);}
73  else if (!strcmp(word, "PM_min" )) {ferr=fscanf(in,"%lf",&physics->Pm );if(ferr!=1)return(1);}
74  else if (!strcmp(word, "lambda" )) {ferr=fscanf(in,"%lf",&physics->l );if(ferr!=1)return(1);}
75  else if (!strcmp(word, "q" )) {ferr=fscanf(in,"%lf",&physics->q );if(ferr!=1)return(1);}
76  else if (!strcmp(word, "tf" )) {ferr=fscanf(in,"%lf",&physics->tf );if(ferr!=1)return(1);}
77  else if (!strcmp(word, "tcl" )) {ferr=fscanf(in,"%lf",&physics->tcl);if(ferr!=1)return(1);}
78  }
79  } else {
80  fprintf(stderr,"Error: Illegal format in file \"physics.inp\".\n");
81  return(1);
82  }
83  }
84  fclose(in);
85 
86  if (!strcmp(solver->SplitHyperbolicFlux,"yes")) {
87  if (!mpi->rank) {
88  fprintf(stderr,"Error in FPPowerSystemInitialize: This physical model does not have a splitting ");
89  fprintf(stderr,"of the hyperbolic term defined.\n");
90  }
91  return(1);
92  }
93 
94  /* initializing physical model-specific functions */
99  solver->Upwind = FPPowerSystemUpwind;
102 
103  /* Calculate and print the PDF integral of the initial solution */
104  IERR FPPowerSystemPostStep(solver->u,solver,mpi,0.0,0); CHECKERR(ierr);
105  if (!mpi->rank) IERR FPPowerSystemPrintStep(solver,mpi,0.0); CHECKERR(ierr);
106 
107  count++;
108  return(0);
109 }
int(* PrintStep)(void *, void *, double)
Definition: hypar.h:344
double FPPowerSystemComputeCFL(void *, void *, double, double)
double FPPowerSystemComputeDiffNumber(void *, void *, double, double)
int(* GFunction)(double *, double *, int, void *, double)
Definition: hypar.h:310
int nvars
Definition: hypar.h:29
#define IERR
Definition: basic.h:16
#define CHECKERR(ierr)
Definition: basic.h:18
char SplitHyperbolicFlux[_MAX_STRING_SIZE_]
Definition: hypar.h:92
double * u
Definition: hypar.h:116
int(* Upwind)(double *, double *, double *, double *, double *, double *, int, void *, double)
Definition: hypar.h:295
int FPPowerSystemPrintStep(void *, void *, double)
int ndims
Definition: hypar.h:26
int(* FFunction)(double *, double *, int, void *, double)
Definition: hypar.h:276
int FPPowerSystemDiffusion(double *, double *, int, void *, double)
double(* ComputeCFL)(void *, void *, double, double)
Definition: hypar.h:269
Structure containing all solver-specific variables and functions.
Definition: hypar.h:23
#define _MAX_STRING_SIZE_
Definition: basic.h:14
#define _MODEL_NDIMS_
Definition: euler1d.h:56
int FPPowerSystemPostStep(double *, void *, void *, double, int)
double(* ComputeDiffNumber)(void *, void *, double, double)
Definition: hypar.h:272
int FPPowerSystemAdvection(double *, double *, int, void *, double)
int(* PostStep)(double *, void *, void *, double, int)
Definition: hypar.h:341
void * physics
Definition: hypar.h:266
Structure of MPI-related variables.
#define _MODEL_NVARS_
Definition: euler1d.h:58
#define _DECLARE_IERR_
Definition: basic.h:17
int FPPowerSystemUpwind(double *, double *, double *, double *, double *, double *, int, void *, double)