HyPar  1.0
Finite-Difference Hyperbolic-Parabolic PDE Solver on Cartesian Grids
Cleanup.c
Go to the documentation of this file.
1 
5 #include <stdio.h>
6 #include <stdlib.h>
7 #include <string.h>
8 #include <basic.h>
9 #include <bandedmatrix.h>
10 #include <tridiagLU.h>
11 #include <boundaryconditions.h>
12 #include <immersedboundaries.h>
13 #include <timeintegration.h>
14 #include <interpolation.h>
15 #include <mpivars.h>
16 #include <simulation_object.h>
17 
18 /* include header files for each physical model */
24 #include <physicalmodels/euler1d.h>
25 #include <physicalmodels/euler2d.h>
28 #include <physicalmodels/numa2d.h>
29 #include <physicalmodels/numa3d.h>
32 #include <physicalmodels/vlasov.h>
33 
34 #if defined(HAVE_CUDA)
35 #include <arrayfunctions_gpu.h>
36 #endif
37 
39 int Cleanup( void *s,
40  int nsims
41  )
42 {
44  int ns;
46 
47  if (nsims == 0) return 0;
48 
49  if (!sim[0].mpi.rank) {
50  printf("Deallocating arrays.\n");
51  }
52 
53  for (ns = 0; ns < nsims; ns++) {
54 
55  if (sim[ns].is_barebones == 1) {
56  fprintf(stderr, "Error in Cleanup(): object is barebones type.\n");
57  return 1;
58  }
59 
60  HyPar* solver = &(sim[ns].solver);
61  MPIVariables* mpi = &(sim[ns].mpi);
62  DomainBoundary* boundary = (DomainBoundary*) solver->boundary;
63  int i;
64 
65  /* Clean up boundary zones */
66  for (i = 0; i < solver->nBoundaryZones; i++) {
67 #if defined(HAVE_CUDA)
68  BCCleanup(&boundary[i], solver->use_gpu);
69 #else
70  BCCleanup(&boundary[i], 0);
71 #endif
72  }
73  free(solver->boundary);
74 
75  /* Clean up immersed boundaries */
76  if (solver->flag_ib) {
77  IERR IBCleanup(solver->ib);
78  free(solver->ib);
79  }
80 
81  /* Clean up any allocations in physical model */
82  if (!strcmp(solver->model,_LINEAR_ADVECTION_DIFFUSION_REACTION_)) {
83  IERR LinearADRCleanup(solver->physics); CHECKERR(ierr);
84  } else if (!strcmp(solver->model,_FP_DOUBLE_WELL_)) {
85  IERR FPDoubleWellCleanup(solver->physics); CHECKERR(ierr);
86  } else if (!strcmp(solver->model,_FP_POWER_SYSTEM_)) {
87  IERR FPPowerSystemCleanup(solver->physics); CHECKERR(ierr);
88  } else if (!strcmp(solver->model,_FP_POWER_SYSTEM_1BUS_)) {
90  } else if (!strcmp(solver->model,_FP_POWER_SYSTEM_3BUS_)) {
92  } else if (!strcmp(solver->model,_EULER_1D_)) {
93  IERR Euler1DCleanup(solver->physics); CHECKERR(ierr);
94  } else if (!strcmp(solver->model,_EULER_2D_)) {
95  IERR Euler2DCleanup(solver->physics); CHECKERR(ierr);
96  } else if (!strcmp(solver->model,_NAVIER_STOKES_2D_)) {
98  } else if (!strcmp(solver->model,_NAVIER_STOKES_3D_)) {
100 #if defined(HAVE_CUDA)
101  if (solver->use_gpu) {
103  }
104 #endif
105  } else if (!strcmp(solver->model,_NUMA2D_)) {
106  IERR Numa2DCleanup(solver->physics); CHECKERR(ierr);
107  } else if (!strcmp(solver->model,_NUMA3D_)) {
108  IERR Numa3DCleanup(solver->physics); CHECKERR(ierr);
109  } else if (!strcmp(solver->model,_SHALLOW_WATER_1D_)) {
110  IERR ShallowWater1DCleanup(solver->physics); CHECKERR(ierr);
111  } else if (!strcmp(solver->model,_SHALLOW_WATER_2D_)) {
112  IERR ShallowWater2DCleanup(solver->physics); CHECKERR(ierr);
113  } else if (!strcmp(solver->model,_VLASOV_)) {
114  IERR VlasovCleanup(solver->physics); CHECKERR(ierr);
115  }
116  free(solver->physics);
117 
118  /* Clean up any allocations from time-integration */
119 #ifdef with_petsc
120  if (!solver->use_petscTS) {
121  if (!strcmp(solver->time_scheme,_RK_)) {
122  IERR TimeExplicitRKCleanup(solver->msti); CHECKERR(ierr);
123  free(solver->msti);
124  } else if (!strcmp(solver->time_scheme,_GLM_GEE_)) {
125  IERR TimeGLMGEECleanup(solver->msti); CHECKERR(ierr);
126  free(solver->msti);
127  }
128  }
129 #else
130  if (!strcmp(solver->time_scheme,_RK_)) {
131  IERR TimeExplicitRKCleanup(solver->msti); CHECKERR(ierr);
132  free(solver->msti);
133  } else if (!strcmp(solver->time_scheme,_GLM_GEE_)) {
134  IERR TimeGLMGEECleanup(solver->msti); CHECKERR(ierr);
135  free(solver->msti);
136  }
137 #endif
138 
139  /* Clean up any spatial reconstruction related allocations */
140  if ( (!strcmp(solver->spatial_scheme_hyp,_FIFTH_ORDER_WENO_ ))
141  || (!strcmp(solver->spatial_scheme_hyp,_FIFTH_ORDER_CRWENO_))
142  || (!strcmp(solver->spatial_scheme_hyp,_FIFTH_ORDER_HCWENO_)) ) {
143 #if defined(HAVE_CUDA)
144  IERR WENOCleanup(solver->interp, solver->use_gpu); CHECKERR(ierr);
145 #else
146  IERR WENOCleanup(solver->interp, 0); CHECKERR(ierr);
147 #endif
148  }
149  if (solver->interp) free(solver->interp);
150  if ( (!strcmp(solver->spatial_scheme_hyp,_FIFTH_ORDER_COMPACT_UPWIND_ ))
151  || (!strcmp(solver->spatial_scheme_hyp,_FIFTH_ORDER_CRWENO_ ))
152  || (!strcmp(solver->spatial_scheme_hyp,_FIFTH_ORDER_HCWENO_ )) ) {
153  IERR CompactSchemeCleanup(solver->compact); CHECKERR(ierr);
154  }
155  if (solver->compact) free(solver->compact);
156  if (solver->lusolver) free(solver->lusolver);
157 
158  /* Free the communicators created */
159  IERR MPIFreeCommunicators(solver->ndims,mpi); CHECKERR(ierr);
160 
161  /* These variables are allocated in Initialize.c */
162  free(solver->dim_global);
163  free(solver->dim_global_ex);
164  free(solver->dim_local);
165  free(solver->index);
166  free(solver->u);
167 #ifdef with_petsc
168  if (solver->u0) free(solver->u0);
169  if (solver->uref) free(solver->uref);
170  if (solver->rhsref) free(solver->rhsref);
171  if (solver->rhs) free(solver->rhs);
172 #endif
173 #ifdef with_librom
174  free(solver->u_rom_predicted);
175 #endif
176  free(solver->iblank);
177  free(solver->x);
178  free(solver->dxinv);
179  free(solver->isPeriodic);
180  free(mpi->iproc);
181  free(mpi->ip);
182  free(mpi->is);
183  free(mpi->ie);
184  free(mpi->bcperiodic);
185  free(mpi->sendbuf);
186  free(mpi->recvbuf);
187  free(solver->VolumeIntegral);
188  free(solver->VolumeIntegralInitial);
189  free(solver->TotalBoundaryIntegral);
190  free(solver->ConservationError);
191  free(solver->stride_with_ghosts);
192  free(solver->stride_without_ghosts);
193 
194 #if defined(HAVE_CUDA)
195  if (solver->use_gpu) {
196  gpuFree(solver->hyp);
197  gpuFree(solver->par);
198  gpuFree(solver->source);
199  gpuFree(solver->uC);
200  gpuFree(solver->fluxC);
201  gpuFree(solver->Deriv1);
202  gpuFree(solver->Deriv2);
203  gpuFree(solver->fluxI);
204  gpuFree(solver->uL);
205  gpuFree(solver->uR);
206  gpuFree(solver->fL);
207  gpuFree(solver->fR);
208  gpuFree(solver->StageBoundaryBuffer);
210  gpuFree(solver->StepBoundaryIntegral);
211 
212  gpuFree(solver->gpu_dim_local);
213  gpuFree(solver->gpu_iblank);
214  gpuFree(solver->gpu_x);
215  gpuFree(solver->gpu_dxinv);
216  gpuFree(solver->gpu_u);
217  } else {
218 #endif
219  free(solver->hyp);
220  free(solver->par);
221  free(solver->source);
222  free(solver->uC);
223  free(solver->fluxC);
224  free(solver->Deriv1);
225  free(solver->Deriv2);
226  free(solver->fluxI);
227  free(solver->uL);
228  free(solver->uR);
229  free(solver->fL);
230  free(solver->fR);
231  free(solver->StageBoundaryIntegral);
232  free(solver->StepBoundaryIntegral);
233 #if defined(HAVE_CUDA)
234  }
235 #endif
236 
237  if (solver->filename_index) free(solver->filename_index);
238 
239  }
240 
241  return(0);
242 }
#define _FP_DOUBLE_WELL_
Definition: fpdoublewell.h:18
int WENOCleanup(void *, int)
Definition: WENOCleanup.c:17
double * Deriv2
Definition: hypar.h:151
#define _FP_POWER_SYSTEM_3BUS_
#define _FIFTH_ORDER_WENO_
Definition: interpolation.h:26
#define _LINEAR_ADVECTION_DIFFUSION_REACTION_
Definition: linearadr.h:21
int Numa3DCleanup(void *)
Definition: Numa3DCleanup.c:5
#define _FIFTH_ORDER_HCWENO_
Definition: interpolation.h:30
int LinearADRCleanup(void *)
double * StageBoundaryBuffer
Definition: hypar.h:462
int NavierStokes2DCleanup(void *)
double * recvbuf
Containts the structures and definitions for boundary condition implementation.
1D Shallow Water Equations
#define IERR
Definition: basic.h:16
MPI related function definitions.
char * filename_index
Definition: hypar.h:197
Contains function definitions for common array operations on GPU.
int NavierStokes3DCleanup(void *)
#define CHECKERR(ierr)
Definition: basic.h:18
Structure defining a simulation.
void * boundary
Definition: hypar.h:159
Header file for TridiagLU.
#define _NAVIER_STOKES_2D_
double * par
Definition: hypar.h:122
Data structure and some function declarations for banded block matrices.
int * dim_global_ex
Definition: hypar.h:75
int BCCleanup(void *, int)
Definition: BCCleanup.c:12
int VlasovCleanup(void *)
Definition: VlasovCleanup.c:10
double * gpu_iblank
Definition: hypar.h:456
double * gpu_dxinv
Definition: hypar.h:458
double * u
Definition: hypar.h:116
int Euler2DCleanup(void *)
Definition: Euler2DCleanup.c:4
void * ib
Definition: hypar.h:443
Some basic definitions and macros.
Contains function declarations for time integration.
#define _VLASOV_
Definition: vlasov.h:37
double * ConservationError
Definition: hypar.h:374
Simulation object.
int flag_ib
Definition: hypar.h:441
double * iblank
Definition: hypar.h:436
int ShallowWater1DCleanup(void *)
double * uR
Definition: hypar.h:139
double * x
Definition: hypar.h:107
double * gpu_x
Definition: hypar.h:457
int use_petscTS
Definition: hypar.h:395
int TimeExplicitRKCleanup(void *)
2D Shallow Water Equations
double * u0
Definition: hypar.h:396
int * stride_with_ghosts
Definition: hypar.h:414
double * fluxI
Definition: hypar.h:136
Structure containing the variables and function pointers defining a boundary.
3D Navier Stokes equations (compressible flows)
int FPPowerSystem1BusCleanup(void *)
int Cleanup(void *s, int nsims)
Definition: Cleanup.c:39
int ndims
Definition: hypar.h:26
double * uC
Definition: hypar.h:131
int gpuNavierStokes3DCleanup(void *)
Vlasov Equation.
#define _FP_POWER_SYSTEM_1BUS_
int FPPowerSystemCleanup(void *)
double * rhs
Definition: hypar.h:399
int ShallowWater2DCleanup(void *)
#define _FP_POWER_SYSTEM_
Definition: fppowersystem.h:44
int * gpu_dim_local
Definition: hypar.h:455
3-Bus Power System model
Structure containing all solver-specific variables and functions.
Definition: hypar.h:23
int * stride_without_ghosts
Definition: hypar.h:416
int Numa2DCleanup(void *)
Definition: Numa2DCleanup.c:5
double * sendbuf
char spatial_scheme_hyp[_MAX_STRING_SIZE_]
Definition: hypar.h:84
#define _SHALLOW_WATER_1D_
double * Deriv1
Definition: hypar.h:151
int MPIFreeCommunicators(int, void *)
double * gpu_u
Definition: hypar.h:459
double * uL
Definition: hypar.h:139
int Euler1DCleanup(void *)
char time_scheme[_MAX_STRING_SIZE_]
Definition: hypar.h:78
#define _SHALLOW_WATER_2D_
2D Navier Stokes equations (compressible flows)
#define _EULER_1D_
Definition: euler1d.h:50
int * index
Definition: hypar.h:102
int TimeGLMGEECleanup(void *)
void * lusolver
Definition: hypar.h:368
Linear Advection-Diffusion-Reaction model.
int * dim_local
Definition: hypar.h:37
void * interp
Definition: hypar.h:362
#define _GLM_GEE_
#define _RK_
int CompactSchemeCleanup(void *)
double * rhsref
Definition: hypar.h:398
#define _NUMA3D_
Definition: numa3d.h:26
double * StageBoundaryIntegral
Definition: hypar.h:382
void * physics
Definition: hypar.h:266
char model[_MAX_STRING_SIZE_]
Definition: hypar.h:263
#define _FIFTH_ORDER_CRWENO_
Definition: interpolation.h:28
double * TotalBoundaryIntegral
Definition: hypar.h:386
int IBCleanup(void *)
Definition: IBCleanup.c:11
double * fR
Definition: hypar.h:139
Structures and function definitions for immersed boundaries.
1D Euler Equations (inviscid, compressible flows)
double * fL
Definition: hypar.h:139
#define _NUMA2D_
Definition: numa2d.h:21
#define _EULER_2D_
Definition: euler2d.h:25
Fokker-Planck Model for 1-Bus Power System.
double * StepBoundaryIntegral
Definition: hypar.h:384
#define _FIFTH_ORDER_COMPACT_UPWIND_
Definition: interpolation.h:24
Structure of MPI-related variables.
int FPDoubleWellCleanup(void *)
int use_gpu
Definition: hypar.h:449
Definitions for the functions computing the interpolated value of the primitive at the cell interface...
void * msti
Definition: hypar.h:366
int nBoundaryZones
Definition: hypar.h:157
double * u_rom_predicted
Definition: hypar.h:403
#define _DECLARE_IERR_
Definition: basic.h:17
double * uref
Definition: hypar.h:397
void gpuFree(void *)
int * dim_global
Definition: hypar.h:33
#define _NAVIER_STOKES_3D_
int FPPowerSystem3BusCleanup(void *)
int * isPeriodic
Definition: hypar.h:162
double * VolumeIntegral
Definition: hypar.h:378
void * compact
Definition: hypar.h:364
double * source
Definition: hypar.h:125
double * hyp
Definition: hypar.h:119
double * VolumeIntegralInitial
Definition: hypar.h:380
double * fluxC
Definition: hypar.h:128
double * dxinv
Definition: hypar.h:110