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

Go to the source code of this file.

Functions

double FPDoubleWellComputeCFL (void *, void *, double, double)
 
double FPDoubleWellComputeDiffNumber (void *, void *, double, double)
 
int FPDoubleWellAdvection (double *, double *, int, void *, double)
 
int FPDoubleWellDiffusion (double *, double *, int, void *, double)
 
int FPDoubleWellUpwind (double *, double *, double *, double *, double *, double *, int, void *, double)
 
int FPDoubleWellPostStep (double *, void *, void *, double, int)
 
int FPDoubleWellPrintStep (void *, void *, double)
 
int FPDoubleWellInitialize (void *s, void *m)
 

Function Documentation

◆ FPDoubleWellComputeCFL()

double FPDoubleWellComputeCFL ( void *  ,
void *  ,
double  ,
double   
)

Definition at line 7 of file FPDoubleWellComputeCFL.c.

8 {
9  HyPar *solver = (HyPar*) s;
10  int d, i, v;
11 
12  int ndims = solver->ndims;
13  int nvars = solver->nvars;
14  int ghosts = solver->ghosts;
15  int *dim = solver->dim_local;
16  double *dxinv = solver->dxinv;
17 
18  int offset = 0;
19  double max_cfl = 0;
20  for (d = 0; d < ndims; d++) {
21  for (i = 0; i < dim[d]; i++) {
22  for (v = 0; v < nvars; v++) {
23  double x; _GetCoordinate_(0,i,dim,ghosts,solver->x,x);
24  double local_cfl = absolute(drift(x)) * dt
25  * dxinv[offset+ghosts+i];
26  if (local_cfl > max_cfl) max_cfl = local_cfl;
27  }
28  }
29  offset += dim[d] + 2*ghosts;
30  }
31 
32  return(max_cfl);
33 }
#define absolute(a)
Definition: math_ops.h:32
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 * dim_local
Definition: hypar.h:37
#define _GetCoordinate_(dir, i, dim, ghosts, x, coord)
Definition: basic.h:31
int ghosts
Definition: hypar.h:52
double * dxinv
Definition: hypar.h:110

◆ FPDoubleWellComputeDiffNumber()

double FPDoubleWellComputeDiffNumber ( void *  ,
void *  ,
double  ,
double   
)

Definition at line 5 of file FPDoubleWellComputeDiffNumber.c.

6 {
7  HyPar *solver = (HyPar*) s;
8  FPDoubleWell *params = (FPDoubleWell*) solver->physics;
9  int d, i, v;
10 
11  int ndims = solver->ndims;
12  int nvars = solver->nvars;
13  int ghosts = solver->ghosts;
14  int *dim = solver->dim_local;
15  double *dxinv = solver->dxinv;
16 
17  int offset = 0;
18  double max_diffno = 0;
19  for (d = 0; d < ndims; d++) {
20  for (i = 0; i < dim[d]; i++) {
21  for (v = 0; v < nvars; v++) {
22  double local_diffno = 0.5 * params->q * dt
23  * dxinv[offset+ghosts+i]
24  * dxinv[offset+ghosts+i];
25  if (local_diffno > max_diffno) max_diffno = local_diffno;
26  }
27  }
28  offset += dim[d] + 2*ghosts;
29  }
30 
31  return(max_diffno);
32 }
int nvars
Definition: hypar.h:29
int ndims
Definition: hypar.h:26
Structure containing all solver-specific variables and functions.
Definition: hypar.h:23
int * dim_local
Definition: hypar.h:37
void * physics
Definition: hypar.h:266
int ghosts
Definition: hypar.h:52
double * dxinv
Definition: hypar.h:110

◆ FPDoubleWellAdvection()

int FPDoubleWellAdvection ( double *  ,
double *  ,
int  ,
void *  ,
double   
)

Definition at line 7 of file FPDoubleWellAdvection.c.

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

◆ FPDoubleWellDiffusion()

int FPDoubleWellDiffusion ( double *  ,
double *  ,
int  ,
void *  ,
double   
)

Definition at line 7 of file FPDoubleWellDiffusion.c.

8 {
9  HyPar *solver = (HyPar*) s;
10  FPDoubleWell *param = (FPDoubleWell*) 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] = 0.5*param->q * 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)
void * physics
Definition: hypar.h:266
int ghosts
Definition: hypar.h:52
#define _ArrayCopy1D_(x, y, size)

◆ FPDoubleWellUpwind()

int FPDoubleWellUpwind ( double *  ,
double *  ,
double *  ,
double *  ,
double *  ,
double *  ,
int  ,
void *  ,
double   
)

Definition at line 7 of file FPDoubleWellUpwind.c.

9 {
10  HyPar *solver = (HyPar*) s;
11  int done,v;
12 
13  int ndims = solver->ndims;
14  int nvars = solver->nvars;
15  int ghosts = solver->ghosts;
16  int *dim = solver->dim_local;
17 
18  int index_outer[ndims], index_inter[ndims], bounds_outer[ndims], bounds_inter[ndims];
19  _ArrayCopy1D_(dim,bounds_outer,ndims); bounds_outer[dir] = 1;
20  _ArrayCopy1D_(dim,bounds_inter,ndims); bounds_inter[dir] += 1;
21 
22  done = 0; _ArraySetValue_(index_outer,ndims,0);
23  while (!done) {
24  _ArrayCopy1D_(index_outer,index_inter,ndims);
25  for (index_inter[dir] = 0; index_inter[dir] < bounds_inter[dir]; index_inter[dir]++) {
26  int p; _ArrayIndex1D_(ndims,bounds_inter,index_inter,0,p);
27  double x1; _GetCoordinate_(0,index_inter[0]-1,dim,ghosts,solver->x,x1);
28  double x2; _GetCoordinate_(0,index_inter[0] ,dim,ghosts,solver->x,x2);
29  double x = 0.5 * ( x1 + x2 );
30  for (v = 0; v < nvars; v++)
31  fI[nvars*p+v] = (drift(x) > 0 ? fL[nvars*p+v] : fR[nvars*p+v] );
32  }
33  _ArrayIncrementIndex_(ndims,bounds_outer,index_outer,done);
34  }
35 
36  return(0);
37 }
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
int ghosts
Definition: hypar.h:52
#define _ArrayCopy1D_(x, y, size)

◆ FPDoubleWellPostStep()

int FPDoubleWellPostStep ( double *  ,
void *  ,
void *  ,
double  ,
int   
)

Definition at line 8 of file FPDoubleWellPostStep.c.

9 {
10  HyPar *solver = (HyPar*) s;
11  MPIVariables *mpi = (MPIVariables*) m;
12  FPDoubleWell *params = (FPDoubleWell*) 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 dx = 1.0 / solver->dxinv[index[0]+ghosts];
25  local_sum += (u[p] * dx);
26  _ArrayIncrementIndex_(ndims,dim,index,done);
27  }
28  double local_integral = local_sum;
29  double global_integral = 0;
30  IERR MPISum_double(&global_integral,&local_integral,1,&mpi->world); CHECKERR(ierr);
31  params->pdf_integral = global_integral;
32 
33  return(0);
34 }
#define IERR
Definition: basic.h:16
#define CHECKERR(ierr)
Definition: basic.h:18
double pdf_integral
Definition: fpdoublewell.h:28
int MPISum_double(double *, double *, int, void *)
Definition: MPISum.c:39
int ndims
Definition: hypar.h:26
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)
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

◆ FPDoubleWellPrintStep()

int FPDoubleWellPrintStep ( void *  ,
void *  ,
double   
)

Definition at line 5 of file FPDoubleWellPrintStep.c.

6 {
7  HyPar *solver = (HyPar*) s;
8  FPDoubleWell *params = (FPDoubleWell*) 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: fpdoublewell.h:28
Structure containing all solver-specific variables and functions.
Definition: hypar.h:23
void * physics
Definition: hypar.h:266

◆ FPDoubleWellInitialize()

int FPDoubleWellInitialize ( void *  s,
void *  m 
)

Definition at line 19 of file FPDoubleWellInitialize.c.

20 {
21  HyPar *solver = (HyPar*) s;
22  MPIVariables *mpi = (MPIVariables*) m;
23  FPDoubleWell *physics = (FPDoubleWell*) solver->physics;
24  int ferr = 0;
26 
27  static int count = 0;
28 
29  if (solver->nvars != _MODEL_NVARS_) {
30  fprintf(stderr,"Error in FPDoubleWellInitializeO(): nvars has to be %d.\n",_MODEL_NVARS_);
31  return(1);
32  }
33  if (solver->ndims != _MODEL_NDIMS_) {
34  fprintf(stderr,"Error in FPDoubleWellInitializeO(): ndims has to be %d.\n",_MODEL_NDIMS_);
35  return(1);
36  }
37 
38  /* reading physical model specific inputs - all processes */
39  if (!mpi->rank) {
40  FILE *in;
41  if (!count) printf("Reading physical model inputs from file \"physics.inp\".\n");
42  in = fopen("physics.inp","r");
43  if (!in) {
44  fprintf(stderr,"Error: File \"physics.inp\" not found.\n");
45  return(1);
46  } else {
47  char word[_MAX_STRING_SIZE_];
48  ferr = fscanf(in,"%s",word); if (ferr != 1) return(1);
49  if (!strcmp(word, "begin")){
50  while (strcmp(word, "end")){
51  ferr = fscanf(in,"%s",word); if (ferr != 1) return(1);
52  if (!strcmp(word, "q")) {
53  /* read diffusion coefficient */
54  ferr = fscanf(in,"%lf",&physics->q);
55  if (ferr != 1) return(1);
56  } else if (strcmp(word,"end")) {
57  char useless[_MAX_STRING_SIZE_];
58  ferr = fscanf(in,"%s",useless); if (ferr != 1) return(ferr);
59  printf("Warning: keyword %s in file \"physics.inp\" with value %s not ",word,useless);
60  printf("recognized or extraneous. Ignoring.\n");
61  }
62  }
63  } else {
64  fprintf(stderr,"Error: Illegal format in file \"physics.inp\".\n");
65  return(1);
66  }
67  }
68  fclose(in);
69  }
70 
71 #ifndef serial
72  IERR MPIBroadcast_double(&physics->q,1,0,&mpi->world); CHECKERR(ierr);
73 #endif
74 
75  if (!strcmp(solver->SplitHyperbolicFlux,"yes")) {
76  if (!mpi->rank) {
77  fprintf(stderr,"Error in FPDoubleWellInitialize: This physical model does not have a splitting ");
78  fprintf(stderr,"of the hyperbolic term defined.\n");
79  }
80  return(1);
81  }
82 
83  /* initializing physical model-specific functions */
88  solver->Upwind = FPDoubleWellUpwind;
91 
92  /* Calculate and print the PDF integral of the initial solution */
93  IERR FPDoubleWellPostStep(solver->u,solver,mpi,0.0,0);CHECKERR(ierr);
94  IERR FPDoubleWellPrintStep(solver,mpi,0.0); CHECKERR(ierr);
95 
96  count++;
97  return(0);
98 }
int(* PrintStep)(void *, void *, double)
Definition: hypar.h:344
int(* GFunction)(double *, double *, int, void *, double)
Definition: hypar.h:310
int nvars
Definition: hypar.h:29
double FPDoubleWellComputeCFL(void *, void *, double, double)
#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 ndims
Definition: hypar.h:26
int(* FFunction)(double *, double *, int, void *, double)
Definition: hypar.h:276
int MPIBroadcast_double(double *, int, int, void *)
Definition: MPIBroadcast.c:9
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
double FPDoubleWellComputeDiffNumber(void *, void *, double, double)
#define _MODEL_NDIMS_
Definition: euler1d.h:56
int FPDoubleWellUpwind(double *, double *, double *, double *, double *, double *, int, void *, double)
int FPDoubleWellDiffusion(double *, double *, int, void *, double)
MPI_Comm world
int FPDoubleWellPrintStep(void *, void *, double)
double(* ComputeDiffNumber)(void *, void *, double, double)
Definition: hypar.h:272
int(* PostStep)(double *, void *, void *, double, int)
Definition: hypar.h:341
void * physics
Definition: hypar.h:266
int FPDoubleWellPostStep(double *, void *, void *, double, int)
Structure of MPI-related variables.
#define _MODEL_NVARS_
Definition: euler1d.h:58
#define _DECLARE_IERR_
Definition: basic.h:17
int FPDoubleWellAdvection(double *, double *, int, void *, double)