PICurv 0.1.0
A Parallel Particle-In-Cell Solver for Curvilinear LES
Loading...
Searching...
No Matches
AnalyticalSolution.h
Go to the documentation of this file.
1#ifndef ANALYTICALSOLUTION_H
2#define ANALYTICALSOLUTION_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 <petscdmcomposite.h>
11
12// Include additional headers
13#include "variables.h" // Shared type definitions
14#include "ParticleSwarm.h" // Particle swarm functions
15#include "walkingsearch.h" // Particle location functions
16#include "grid.h" // Grid functions
17#include "logging.h" // Logging macros
18#include "io.h" // Data Input and Output functions
19#include "setup.h" // Setup Module required for array allocation and deallocation
20#include "interpolation.h" // All the different interpolation routines required
21
22/**
23 * @brief Generic macro to select the appropriate local Cartesian field setter based on field type.
24 *
25 * This macro dispatches to:
26 * - SetLocalCartesianField_Vector if the field value pointer is of type Cmpnts*
27 * - SetLocalCartesianField_Scalar if the field value pointer is of type PetscReal*
28 *
29 * @param fieldName Field name for logging.
30 * @param fieldValue Pointer to the cell's field value (either Cmpnts* or PetscReal*).
31 * @param coor Pointer to the coordinate (Cmpnts*) to be used in the computation.
32 */
33#define SetLocalCartesianField(fieldName, fieldValue, coor, FieldInitialization, Constant) _Generic((fieldValue), \
34 Cmpnts*: SetLocalCartesianField_Vector, \
35 PetscReal*: SetLocalCartesianField_Scalar \
36 )(fieldName, fieldValue, coor,FieldInitialization,Constant)
37
38
39
40// function declarations //
41
42
43/**
44 * @brief Sets the local Cartesian scalar field based on input coordinates.
45 *
46 * This function computes the scalar field value by combining the sine of the input
47 * coordinate values. In this example, the scalar field is computed as the sum of the
48 * sine functions of the x, y, and z coordinates.
49 *
50 * @param[in] fieldName Pointer to a string representing the field name (for logging purposes).
51 * @param[in,out] scalarField Pointer to the PetscReal where the computed scalar field value will be stored.
52 * @param[in] coor Pointer to the `Cmpnts` structure containing the input coordinate values.
53 *
54 * @return PetscErrorCode Returns 0 on successful execution, non-zero on failure.
55 */
56PetscErrorCode SetLocalCartesianField_Scalar(const char *fieldName, PetscReal *scalarField, Cmpnts *coor,PetscInt FieldInitialization, PetscReal Constant);
57
58/**
59 * @brief Sets the local Cartesian vector field based on input coordinates.
60 *
61 * This function computes the vector components by applying the sine function to the
62 * input coordinate values along the x, y, and z directions.
63 *
64 * @param[in] fieldName Pointer to a string representing the field name (for logging purposes).
65 * @param[in,out] vecField Pointer to the `Cmpnts` structure where the computed vector field will be stored.
66 * @param[in] coor Pointer to the `Cmpnts` structure containing the input coordinate values.
67 *
68 * @return PetscErrorCode Returns 0 on successful execution, non-zero on failure.
69 */
70PetscErrorCode SetLocalCartesianField_Vector(const char *fieldName, Cmpnts *vecField, Cmpnts *coor, PetscInt FieldInitialization, PetscReal Constant);
71
72/**
73 * @brief Sets an analytical Cartesian field (scalar or vector) for cell centers based on a field name.
74 *
75 * This function looks up the field within the user context by comparing the provided field name
76 * against the supported names:
77 * - Vector fields: "Ucat" and "Ucont"
78 * - Scalar fields: "P" and "nvert"
79 *
80 * If the field is found, the function verifies that the DM block size matches the expected value
81 * (3 for vector fields and 1 for scalar fields), retrieves local DM information for both the coordinate DM (da)
82 * and the cell-centered DM (fda), interpolates the corner-based coordinates (from da) to cell centers,
83 * and then updates the field using the generic helper macro SetLocalCartesianField. Interior cells are
84 * updated with the interpolated coordinates, and boundary cells are updated using the original coordinate data.
85 *
86 * If the field name is not found in the user context, the function throws an error.
87 *
88 * @param[in] user Pointer to the UserCtx structure containing:
89 * - da: DM for coordinate (corner) data.
90 * - fda: DM for cell-centered data.
91 * - Ucat: vector field (Cmpnts ***)
92 * - Ucont: vector field (Cmpnts ***)
93 * - P: scalar field (PetscReal ***)
94 * - nvert: scalar field (PetscReal ***)
95 * @param[in] fieldName Name of the field to update.
96 *
97 * @return PetscErrorCode 0 on success, non-zero on failure.
98 */
99PetscErrorCode SetAnalyticalCartesianField(UserCtx *user, const char *fieldName);
100
101/**
102 * @brief Sets an analytical contravariant field (Ucont) at nodes, where each component
103 * is evaluated using the physical coordinates of its respective face center.
104 *
105 * This function:
106 * 1. Interpolates nodal physical coordinates to the centers of i-faces, j-faces, and k-faces
107 * for the cells relevant to the current rank, storing them in temporary local arrays.
108 * 2. For each node (i,j,k) owned by the current rank (where Ucont components are stored):
109 * a. Retrieves the pre-calculated physical coordinate of the i-face center.
110 * b. Evaluates the analytical function for U^1 using this i-face center coordinate.
111 * c. Retrieves the pre-calculated physical coordinate of the j-face center.
112 * d. Evaluates the analytical function for U^2 using this j-face center coordinate.
113 * e. Retrieves the pre-calculated physical coordinate of the k-face center.
114 * f. Evaluates the analytical function for U^3 using this k-face center coordinate.
115 * g. Stores these (U^1, U^2, U^3) into user->Ucont[k][j][i].
116 *
117 * @param[in] user Pointer to the UserCtx structure.
118 * @param[in] fieldName Name of the field to update (currently only "Ucont").
119 *
120 * @return PetscErrorCode 0 on success, non-zero on failure.
121 */
122PetscErrorCode SetAnalyticalContravariantField(UserCtx *user, const char *fieldName);
123
124
125/**
126 * @brief Applies the analytical solution to the position vector.
127 *
128 * This function updates each entry in the provided PETSc vector by computing its sine,
129 * thereby replacing each position with sin(position).
130 *
131 * @param tempVec The PETSc Vec containing particle positions which will be used to store the velocities.
132 * @return PetscErrorCode Returns 0 on success.
133 */
134PetscErrorCode SetAnalyticalSolution(Vec tempVec, PetscInt FieldInitialization);
135
136/**
137 * @brief Applies analytical boundary conditions to a specified global vector or scalar field.
138 *
139 * This function acts as a dispatcher based on the fieldName. It determines the
140 * appropriate DM and Vec from the UserCtx and calls field-specific helper routines
141 * (e.g., ApplyAnalyticalBC_Vector, ApplyAnalyticalBC_Scalar) to set boundary values
142 * according to a predefined analytical function (e.g., sin(coord)).
143 *
144 * This should be called AFTER setting interior values and BEFORE updating local ghosts.
145 *
146 * @param[in] user Pointer to the UserCtx containing DMs, coordinates, and field Vecs.
147 * @param[in] fieldName Name of the field to apply BCs to ("Ucat", "Ucont", "P", "Nvert").
148 *
149 * @return PetscErrorCode 0 on success, non-zero on failure.
150 */
151PetscErrorCode ApplyAnalyticalBC(UserCtx *user, const char *fieldName);
152
153/**
154 * @brief Sets the initial values for the INTERIOR of a specified Eulerian field.
155 *
156 * This function initializes the interior nodes of `Ucont` based on a profile selected
157 * by `user->FieldInitialization`. It explicitly skips any node that lies on a global
158 * boundary, as those values are set by the Boundary System's `Initialize` methods.
159 *
160 * The initialization is directional, aligned with the primary INLET face that was
161 * identified by the parser. This ensures the initial flow is physically meaningful.
162 *
163 * Supported `user->FieldInitialization` profiles for "Ucont":
164 * - 0: Zero Velocity. All interior components of Ucont are set to 0.
165 * - 1: Constant Normal Velocity. The contravariant velocity component normal to the
166 * inlet direction is set such that the physical velocity normal to those grid
167 * planes is a constant `uin`. Other contravariant components are zero.
168 * - 2: Poiseuille Normal Velocity. The contravariant component normal to the
169 * inlet direction is set with a parabolic profile.
170 *
171 * @param user The main UserCtx struct, containing all simulation data and configuration.
172 * @param fieldName A string ("Ucont" or "P") identifying which field to initialize.
173 * @return PetscErrorCode 0 on success.
174 */
175PetscErrorCode SetInitialInteriorField(UserCtx *user, const char *fieldName);
176
177#endif // ANALYTICALSOLUTION_H
PetscErrorCode SetAnalyticalCartesianField(UserCtx *user, const char *fieldName)
Sets an analytical Cartesian field (scalar or vector) for cell centers based on a field name.
PetscErrorCode SetAnalyticalContravariantField(UserCtx *user, const char *fieldName)
Sets an analytical contravariant field (Ucont) at nodes, where each component is evaluated using the ...
PetscErrorCode SetInitialInteriorField(UserCtx *user, const char *fieldName)
Sets the initial values for the INTERIOR of a specified Eulerian field.
PetscErrorCode SetAnalyticalSolution(Vec tempVec, PetscInt FieldInitialization)
Applies the analytical solution to the position vector.
PetscErrorCode SetLocalCartesianField_Scalar(const char *fieldName, PetscReal *scalarField, Cmpnts *coor, PetscInt FieldInitialization, PetscReal Constant)
Sets the local Cartesian scalar field based on input coordinates.
PetscErrorCode ApplyAnalyticalBC(UserCtx *user, const char *fieldName)
Applies analytical boundary conditions to a specified global vector or scalar field.
PetscErrorCode SetLocalCartesianField_Vector(const char *fieldName, Cmpnts *vecField, Cmpnts *coor, PetscInt FieldInitialization, PetscReal Constant)
Sets the local Cartesian vector field based on input coordinates.
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.
A 3D point or vector with PetscScalar components.
Definition variables.h:99
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.