PICurv
0.1.0
A Parallel Particle-In-Cell Solver for Curvilinear LES
Loading...
Searching...
No Matches
initialcondition.h
Go to the documentation of this file.
1
#ifndef INITIALCONDITION_H
2
#define INITIALCONDITION_H
3
4
#include <petscpf.h>
5
#include <petscdmswarm.h>
6
#include <stdlib.h>
7
#include <time.h>
8
#include <math.h>
9
#include <petsctime.h>
10
#include <petscsys.h>
11
#include <petscdmcomposite.h>
12
#include <petscsystypes.h>
13
// Include additional headers
14
#include "
variables.h
"
// Shared type definitions
15
#include "
ParticleSwarm.h
"
// Particle swarm functions
16
#include "
walkingsearch.h
"
// Particle location functions
17
#include "
grid.h
"
// Grid functions
18
#include "
logging.h
"
// Logging macros
19
#include "
io.h
"
// Data Input and Output functions
20
#include "
interpolation.h
"
// Interpolation routines
21
#include "
AnalyticalSolution.h
"
// Analytical Solution for testing
22
#include "
ParticleMotion.h
"
// Functions related to motion of particles
23
#include "
Boundaries.h
"
// Functions related to Boundary condition
24
#include "
simulation.h
"
25
26
/**
27
* @brief Sets the initial values for the INTERIOR of a specified Eulerian field.
28
*
29
* This function initializes the interior nodes of `Ucont` based on a profile selected
30
* by `user->FieldInitialization`. It explicitly skips any node that lies on a global
31
* boundary, as those values are set by the Boundary System's `Initialize` methods.
32
*
33
* The initialization is directional, aligned with the primary INLET face that was
34
* identified by the parser. This ensures the initial flow is physically meaningful.
35
*
36
* Supported `user->FieldInitialization` profiles for "Ucont":
37
* - 0: Zero Velocity. All interior components of Ucont are set to 0.
38
* - 1: Constant Normal Velocity. The contravariant velocity component normal to the
39
* inlet direction is set such that the physical velocity normal to those grid
40
* planes is a constant `uin`. Other contravariant components are zero.
41
* - 2: Poiseuille Normal Velocity. The contravariant component normal to the
42
* inlet direction is set with a parabolic profile.
43
*
44
* @param user The main UserCtx struct, containing all simulation data and configuration.
45
* @param fieldName A string ("Ucont" or "P") identifying which field to initialize.
46
* @return PetscErrorCode 0 on success.
47
*/
48
PetscErrorCode
SetInitialInteriorField
(
UserCtx
*user,
const
char
*fieldName);
49
50
/**
51
* @brief High-level orchestrator to set the complete initial state of the Eulerian solver.
52
*
53
* This function is called once from main() before the time loop begins. It inspects
54
* the simulation context to determine whether to perform a fresh start (t=0) or
55
* restart from saved files. It then delegates to the appropriate helper function.
56
* Finally, it initializes the solver's history vectors (Ucont_o, P_o, etc.)
57
* to ensure the first time step has the necessary data.
58
*/
59
PetscErrorCode
InitializeEulerianState
(
SimCtx
*simCtx);
60
61
#endif
// INITIALCONDITION_H
AnalyticalSolution.h
Boundaries.h
ParticleMotion.h
Header file for Particle Motion and migration related functions.
ParticleSwarm.h
Header file for Particle Swarm management functions.
grid.h
Public interface for grid, solver, and metric setup routines.
SetInitialInteriorField
PetscErrorCode SetInitialInteriorField(UserCtx *user, const char *fieldName)
Sets the initial values for the INTERIOR of a specified Eulerian field.
Definition
initialcondition.c:32
InitializeEulerianState
PetscErrorCode InitializeEulerianState(SimCtx *simCtx)
High-level orchestrator to set the complete initial state of the Eulerian solver.
Definition
initialcondition.c:381
interpolation.h
io.h
Public interface for data input/output routines.
logging.h
Logging utilities and macros for PETSc-based applications.
simulation.h
variables.h
Main header file for a complex fluid dynamics solver.
SimCtx
The master context for the entire simulation.
Definition
variables.h:538
UserCtx
User-defined context containing data specific to a single computational grid level.
Definition
variables.h:661
walkingsearch.h
Header file for particle location functions using the walking search algorithm.
include
initialcondition.h
Generated by
1.9.8