PICurv 0.1.0
A Parallel Particle-In-Cell Solver for Curvilinear LES
Loading...
Searching...
No Matches
Boundaries.h
Go to the documentation of this file.
1#ifndef BOUNDARIES_H
2#define BOUNDARIES_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
14// Include additional headers
15#include "variables.h" // Shared type definitions
16#include "ParticleSwarm.h" // Particle swarm functions
17#include "walkingsearch.h" // Particle location functions
18#include "grid.h" // Grid functions
19#include "logging.h" // Logging macros
20#include "io.h" // Data Input and Output functions
21#include "interpolation.h" // Interpolation routines
22#include "AnalyticalSolution.h" // Analytical Solution for testing
23#include "ParticleMotion.h" // Functions related to motion of particles
24#include "BC_Handlers.h" // Boundary Handlers
25#include "wallfunction.h" // wall functions for LES
26//================================================================================
27//
28// PUBLIC SYSTEM-LEVEL FUNCTIONS
29//
30// These are the main entry points for interacting with the boundary system.
31//
32//================================================================================
33
34/**
35 * @brief Initializes the entire boundary system.
36 * @param user The main UserCtx struct.
37 * @param bcs_filename The path to the boundary conditions configuration file.
38 */
39PetscErrorCode BoundarySystem_Create(UserCtx *user, const char *bcs_filename);
40
41/**
42 * @brief Executes one full boundary condition update cycle for a time step.
43 * @param user The main UserCtx struct.
44 */
45PetscErrorCode BoundarySystem_ExecuteStep(UserCtx *user);
46
47/**
48 * @brief Cleans up and destroys all boundary system resources.
49 * @param user The main UserCtx struct.
50 */
51PetscErrorCode BoundarySystem_Destroy(UserCtx *user);
52
53/**
54 * @brief Determines if the current MPI rank owns any part of the globally defined inlet face,
55 * making it responsible for placing particles on that portion of the surface.
56 *
57 * The determination is based on the rank's owned nodes (from `DMDALocalInfo`) and
58 * the global node counts, in conjunction with the `user->identifiedInletBCFace`.
59 * A rank can service an inlet face if it owns the cells adjacent to that global boundary
60 * and has a non-zero extent (owns cells) in the tangential dimensions of that face.
61 *
62 * @param user Pointer to the UserCtx structure, containing `identifiedInletBCFace`.
63 * @param info Pointer to the DMDALocalInfo for the current rank's DA (node-based).
64 * @param IM_nodes_global Global number of nodes in the I-direction (e.g., user->IM + 1 if user->IM is cell count).
65 * @param JM_nodes_global Global number of nodes in the J-direction.
66 * @param KM_nodes_global Global number of nodes in the K-direction.
67 * @param[out] can_service_inlet_out Pointer to a PetscBool; set to PETSC_TRUE if the rank
68 * services (part of) the inlet, PETSC_FALSE otherwise.
69 * @return PetscErrorCode 0 on success, non-zero on failure.
70 */
71PetscErrorCode CanRankServiceInletFace(UserCtx *user, const DMDALocalInfo *info,
72 PetscInt IM_nodes_global, PetscInt JM_nodes_global, PetscInt KM_nodes_global,
73 PetscBool *can_service_inlet_out);
74
75/**
76 * @brief Determines if the current MPI rank owns any part of a specified global face.
77 *
78 * This function is a general utility for parallel boundary operations. It checks if the
79 * local domain of the current MPI rank is adjacent to a specified global boundary face.
80 * A rank "services" a face if it owns the cells adjacent to that global boundary and has
81 * a non-zero extent (i.e., owns at least one cell) in the tangential dimensions of that face.
82 *
83 * @param info Pointer to the DMDALocalInfo for the current rank's DA.
84 * @param face_id The specific global face (e.g., BC_FACE_NEG_Z) to check.
85 * @param[out] can_service_out Pointer to a PetscBool; set to PETSC_TRUE if the rank
86 * services the face, PETSC_FALSE otherwise.
87 * @return PetscErrorCode 0 on success.
88 */
89PetscErrorCode CanRankServiceFace(const DMDALocalInfo *info, BCFace face_id, PetscBool *can_service_out);
90
91
92/**
93 * @brief Assuming the current rank services the inlet face, this function selects a random
94 * cell (owned by this rank on that face) and random logical coordinates within that cell,
95 * suitable for placing a particle on the inlet surface.
96 *
97 * It is the caller's responsibility to ensure CanRankServiceInletFace returned true.
98 *
99 * @param user Pointer to UserCtx.
100 * @param info Pointer to DMDALocalInfo for the current rank (node-based).
101 * @param xs_gnode, ys_gnode, zs_gnode Local starting node indices (incl. ghosts) for the rank's DA.
102 * @param IM_nodes_global, JM_nodes_global, KM_nodes_global Global node counts.
103 * @param rand_logic_i_ptr, rand_logic_j_ptr, rand_logic_k_ptr Pointers to RNGs for logical coords.
104 * @param[out] ci_metric_lnode_out, cj_metric_lnode_out, ck_metric_lnode_out Local node indices of the selected cell's origin (these are local to the rank's DA including ghosts).
105 * @param[out] xi_metric_logic_out, eta_metric_logic_out, zta_metric_logic_out Logical coords [0,1] within the cell.
106 * @return PetscErrorCode
107 */
109 UserCtx *user, const DMDALocalInfo *info,
110 PetscInt xs_gnode_rank, PetscInt ys_gnode_rank, PetscInt zs_gnode_rank, // Local starting node index (with ghosts) of the rank's DA patch
111 PetscInt IM_nodes_global, PetscInt JM_nodes_global, PetscInt KM_nodes_global,
112 PetscRandom *rand_logic_i_ptr, PetscRandom *rand_logic_j_ptr, PetscRandom *rand_logic_k_ptr,
113 PetscInt *ci_metric_lnode_out, PetscInt *cj_metric_lnode_out, PetscInt *ck_metric_lnode_out,
114 PetscReal *xi_metric_logic_out, PetscReal *eta_metric_logic_out, PetscReal *zta_metric_logic_out);
115
116PetscErrorCode TranslateModernBCsToLegacy(UserCtx *user);
117
118/* Boundary condition defination (array user->bctype[0-5]):
119 0: interpolation/interface
120 -1: wallfunction
121 1: solid wall (not moving)
122 2: moving solid wall (U=1)
123 3: slip wall/symmetry
124 5: Inlet
125 4: Outlet
126 6: farfield
127 7: periodic
128 8: Characteristic BC
129 9: Analytical Vortex
130 10: Oulet Junction
131 11: Annulus
132 12: Ogrid
133 13: Rheology
134 14: Outlet with Interface
135
136*/
137
138PetscErrorCode InflowFlux(UserCtx *user);
139
140PetscErrorCode OutflowFlux(UserCtx *user);
141
142PetscErrorCode FormBCS(UserCtx *user);
143
144/**
145 * @brief Acts as a temporary bridge to the legacy boundary condition implementation.
146 *
147 * This function is the key to our integration strategy. It matches the signature
148 * of the modern `BoundarySystem_ExecuteStep` function that SetEulerianFields
149 * expects to call.
150 *
151 * However, instead of containing new handler-based logic, it simply calls the
152 * monolithic legacy `FormBCS` function. This allows the modern orchestrator to
153 * drive the old solver logic without modification.
154 *
155 * @param user The UserCtx for a single block.
156 * @return PetscErrorCode
157 */
158PetscErrorCode BoundarySystem_ExecuteStep_Legacy(UserCtx *user);
159
160#endif // BOUNDARIES_H
PetscErrorCode BoundarySystem_ExecuteStep_Legacy(UserCtx *user)
Acts as a temporary bridge to the legacy boundary condition implementation.
Definition Boundaries.c:440
PetscErrorCode CanRankServiceFace(const DMDALocalInfo *info, BCFace face_id, PetscBool *can_service_out)
Determines if the current MPI rank owns any part of a specified global face.
Definition Boundaries.c:130
PetscErrorCode OutflowFlux(UserCtx *user)
PetscErrorCode BoundarySystem_Create(UserCtx *user, const char *bcs_filename)
Initializes the entire boundary system.
PetscErrorCode CanRankServiceInletFace(UserCtx *user, const DMDALocalInfo *info, PetscInt IM_nodes_global, PetscInt JM_nodes_global, PetscInt KM_nodes_global, PetscBool *can_service_inlet_out)
Determines if the current MPI rank owns any part of the globally defined inlet face,...
Definition Boundaries.c:24
PetscErrorCode InflowFlux(UserCtx *user)
Definition Boundaries.c:819
PetscErrorCode GetRandomCellAndLogicOnInletFace(UserCtx *user, const DMDALocalInfo *info, PetscInt xs_gnode_rank, PetscInt ys_gnode_rank, PetscInt zs_gnode_rank, PetscInt IM_nodes_global, PetscInt JM_nodes_global, PetscInt KM_nodes_global, PetscRandom *rand_logic_i_ptr, PetscRandom *rand_logic_j_ptr, PetscRandom *rand_logic_k_ptr, PetscInt *ci_metric_lnode_out, PetscInt *cj_metric_lnode_out, PetscInt *ck_metric_lnode_out, PetscReal *xi_metric_logic_out, PetscReal *eta_metric_logic_out, PetscReal *zta_metric_logic_out)
Assuming the current rank services the inlet face, this function selects a random cell (owned by this...
Definition Boundaries.c:225
PetscErrorCode BoundarySystem_ExecuteStep(UserCtx *user)
Executes one full boundary condition update cycle for a time step.
PetscErrorCode BoundarySystem_Destroy(UserCtx *user)
Cleans up and destroys all boundary system resources.
PetscErrorCode FormBCS(UserCtx *user)
PetscErrorCode TranslateModernBCsToLegacy(UserCtx *user)
Definition Boundaries.c:391
Header file for Particle Motion and migration related functions.
Public interface for grid, solver, and metric setup routines.
Public interface for data input/output routines.
Logging utilities and macros for PETSc-based applications.
Main header file for a complex fluid dynamics solver.
BCFace
Identifies the six logical faces of a structured computational block.
Definition variables.h:199
User-defined context containing data specific to a single computational grid level.
Definition variables.h:630
Header file for particle location functions using the walking search algorithm.