PICurv 0.1.0
A Parallel Particle-In-Cell Solver for Curvilinear LES
Loading...
Searching...
No Matches
initialcondition.h File Reference
#include <petscpf.h>
#include <petscdmswarm.h>
#include <stdlib.h>
#include <time.h>
#include <math.h>
#include <petsctime.h>
#include <petscsys.h>
#include <petscdmcomposite.h>
#include <petscsystypes.h>
#include "variables.h"
#include "ParticleSwarm.h"
#include "walkingsearch.h"
#include "grid.h"
#include "logging.h"
#include "io.h"
#include "interpolation.h"
#include "AnalyticalSolution.h"
#include "ParticleMotion.h"
#include "Boundaries.h"
#include "simulation.h"
Include dependency graph for initialcondition.h:
This graph shows which files directly or indirectly include this file:

Go to the source code of this file.

Functions

PetscErrorCode SetInitialInteriorField (UserCtx *user, const char *fieldName)
 Sets the initial values for the INTERIOR of a specified Eulerian field.
 
PetscErrorCode InitializeEulerianState (SimCtx *simCtx)
 High-level orchestrator to set the complete initial state of the Eulerian solver.
 

Function Documentation

◆ SetInitialInteriorField()

PetscErrorCode SetInitialInteriorField ( UserCtx user,
const char *  fieldName 
)

Sets the initial values for the INTERIOR of a specified Eulerian field.

This function initializes the interior nodes of Ucont based on a profile selected by user->FieldInitialization. It explicitly skips any node that lies on a global boundary, as those values are set by the Boundary System's Initialize methods.

The initialization is directional, aligned with the primary INLET face that was identified by the parser. This ensures the initial flow is physically meaningful.

Supported user->FieldInitialization profiles for "Ucont":

  • 0: Zero Velocity. All interior components of Ucont are set to 0.
  • 1: Constant Normal Velocity. The contravariant velocity component normal to the inlet direction is set such that the physical velocity normal to those grid planes is a constant uin. Other contravariant components are zero.
  • 2: Poiseuille Normal Velocity. The contravariant component normal to the inlet direction is set with a parabolic profile.
Parameters
userThe main UserCtx struct, containing all simulation data and configuration.
fieldNameA string ("Ucont" or "P") identifying which field to initialize.
Returns
PetscErrorCode 0 on success.

Definition at line 30 of file initialcondition.c.

31{
32 PetscErrorCode ierr;
33 PetscFunctionBeginUser;
34
35 SimCtx *simCtx = user->simCtx;
36
37 LOG_ALLOW(GLOBAL, LOG_INFO, "Setting initial INTERIOR field for '%s' with profile %d.\n", fieldName, simCtx->FieldInitialization);
38
39 // This function currently only implements logic for Ucont.
40 if (strcmp(fieldName, "Ucont") != 0) {
41 LOG_ALLOW(GLOBAL, LOG_DEBUG, "Skipping SetInitialInteriorField for non-Ucont field '%s'.\n", fieldName);
42 PetscFunctionReturn(0);
43 }
44
45 // --- 1. Get references to required data and PETSc arrays ---
46 DM fieldDM = user->fda;
47 Vec fieldVec = user->Ucont;
48 DMDALocalInfo info;
49 ierr = DMDAGetLocalInfo(fieldDM, &info); CHKERRQ(ierr);
50
51 Vec localCoor;
52 Cmpnts ***coor_arr;
53 ierr = DMGetCoordinatesLocal(user->da, &localCoor); CHKERRQ(ierr);
54 ierr = DMDAVecGetArrayRead(user->fda, localCoor, &coor_arr); CHKERRQ(ierr);
55
56 Cmpnts ***csi_arr, ***eta_arr, ***zet_arr;
57 ierr = DMDAVecGetArrayRead(user->fda, user->lCsi, &csi_arr); CHKERRQ(ierr);
58 ierr = DMDAVecGetArrayRead(user->fda, user->lEta, &eta_arr); CHKERRQ(ierr);
59 ierr = DMDAVecGetArrayRead(user->fda, user->lZet, &zet_arr); CHKERRQ(ierr);
60
61 // --- 2. Compute Cell-Center Coordinates (only if needed by the selected profile) ---
62 Cmpnts ***cent_coor = NULL;
63 PetscInt xs_cell=0, xm_cell=0, ys_cell=0, ym_cell=0, zs_cell=0, zm_cell=0;
64
65 if (simCtx->FieldInitialization == 2) { // Profile 2 (Poiseuille) requires cell centers.
66 /*
67 ierr = GetOwnedCellRange(&info, 0, &xs_cell, &xm_cell); CHKERRQ(ierr);
68 ierr = GetOwnedCellRange(&info, 1, &ys_cell, &ym_cell); CHKERRQ(ierr);
69 ierr = GetOwnedCellRange(&info, 2, &zs_cell, &zm_cell); CHKERRQ(ierr);
70 */
71 ierr = GetGhostedCellRange(&info,&user->neighbors,0, &xs_cell, &xm_cell); CHKERRQ(ierr);
72 ierr = GetGhostedCellRange(&info,&user->neighbors,1, &ys_cell, &ym_cell); CHKERRQ(ierr);
73 ierr = GetGhostedCellRange(&info,&user->neighbors,2, &zs_cell, &zm_cell); CHKERRQ(ierr);
74
75 if (xm_cell > 0 && ym_cell > 0 && zm_cell > 0) {
76 ierr = Allocate3DArray(&cent_coor, zm_cell, ym_cell, xm_cell); CHKERRQ(ierr);
77 ierr = InterpolateFieldFromCornerToCenter(coor_arr, cent_coor, user); CHKERRQ(ierr);
78 LOG_ALLOW(LOCAL, LOG_DEBUG, "Computed temporary cell-center coordinates for Poiseuille profile.\n");
79 }
80 }
81
82 // --- 3. Loop Over Owned Nodes and Apply Initial Condition to Interior ---
83 Cmpnts ***ucont_arr;
84 ierr = DMDAVecGetArray(fieldDM, fieldVec, &ucont_arr); CHKERRQ(ierr);
85
86 PetscInt i, j, k;
87 const PetscInt mx = info.mx, my = info.my, mz = info.mz; // Global node dimensions
88 const PetscInt xs = info.xs, xe = info.xs + info.xm;
89 const PetscInt ys = info.ys, ye = info.ys + info.ym;
90 const PetscInt zs = info.zs, ze = info.zs + info.zm;
91
92 const Cmpnts uin = simCtx->InitialConstantContra; // Max/average velocity from user options
93 // Flow into a negative face (e.g., -Zeta at k=0) is in the positive physical direction (+z).
94 // Flow into a positive face (e.g., +Zeta at k=mz-1) is in the negative physical direction (-z).
95
96 LOG_ALLOW(GLOBAL,LOG_DEBUG,"InitialConstantContra = (%.3f, %.3f, %.3f)\n",(double)uin.x, (double)uin.y, (double)uin.z);
97
98 const PetscReal flow_direction_sign = (user->identifiedInletBCFace % 2 == 0) ? 1.0 : -1.0;
99
100 for (k = zs; k < ze; k++) {
101 for (j = ys; j < ye; j++) {
102 for (i = xs; i < xe; i++) {
103
104 // The crucial check to ensure we only modify interior nodes.
105 const PetscBool is_interior = (i > 0 && i < mx - 1 &&
106 j > 0 && j < my - 1 &&
107 k > 0 && k < mz - 1);
108
109 if (is_interior) {
110 Cmpnts ucont_val = {0.0, 0.0, 0.0}; // Default to zero velocity
111 PetscReal normal_velocity_mag = 0.0;
112
113 // Step A: Determine the magnitude of the desired physical normal velocity.
114 switch (simCtx->FieldInitialization) {
115 case 0: // Zero initial velocity
116 normal_velocity_mag = 0.0;
117 break;
118 case 1: /* Constant Normal Velocity */
121 normal_velocity_mag = uin.x;
123 normal_velocity_mag = uin.y;
124 } else {
125 normal_velocity_mag = uin.z;
126 }
127 break;
128 case 2: // Poiseuille Normal Velocity
129 {
130 // This profile assumes flow is aligned with the k-index direction.
131 // It uses grid indices (i,j) to define the cross-section, which works for bent geometries.
132 PetscReal u0 = 0.0;
133 if (user->identifiedInletBCFace <= BC_FACE_POS_X) u0 = uin.x;
134 else if (user->identifiedInletBCFace <= BC_FACE_POS_Y) u0 = uin.y;
135 else u0 = uin.z; // Assumes Z-like inlet direction
136
137 // Define channel geometry in "index space" based on global grid dimensions
138 // We subtract 2.0 because the interior runs from index 1 to mx-2 (or my-2).
139 const PetscReal i_width = (PetscReal)(mx - 2);
140 const PetscReal j_width = (PetscReal)(my - 2);
141 const PetscReal i_center = 1.0 + i_width / 2.0;
142 const PetscReal j_center = 1.0 + j_width / 2.0;
143
144 // Create normalized coordinates for the current point (i,j), ranging from -1 to 1
145 const PetscReal i_norm = (i - i_center) / (i_width / 2.0);
146 const PetscReal j_norm = (j - j_center) / (j_width / 2.0);
147
148 // Apply the parabolic profile for a rectangular/square channel
149 // V(i,j) = V_max * (1 - i_norm^2) * (1 - j_norm^2)
150 const PetscReal profile_i = 1.0 - i_norm * i_norm;
151 const PetscReal profile_j = 1.0 - j_norm * j_norm;
152 normal_velocity_mag = u0 * profile_i * profile_j;
153
154 // Clamp to zero for any points outside the channel (or due to minor float errors)
155 if (normal_velocity_mag < 0.0) {
156 normal_velocity_mag = 0.0;
157 }
158 }
159 break;
160 /*
161 PetscReal r_sq = 0.0;
162 const PetscInt i_local = i - xs_cell, j_local = j - ys_cell, k_local = k - zs_cell;
163 if (cent_coor && i_local >= 0 && i_local < xm_cell && j_local >= 0 && j_local < ym_cell && k_local >= 0 && k_local < zm_cell) {
164 const Cmpnts* center = &cent_coor[k_local][j_local][i_local];
165 if (user->identifiedInletBCFace <= BC_FACE_POS_X) r_sq = center->y * center->y + center->z * center->z;
166 else if (user->identifiedInletBCFace <= BC_FACE_POS_Y) r_sq = center->x * center->x + center->z * center->z;
167 else r_sq = center->x * center->x + center->y * center->y;
168 // pick the correct contravariant component for center‐line speed
169 PetscReal u0;
170 if (user->identifiedInletBCFace == BC_FACE_NEG_X ||
171 user->identifiedInletBCFace == BC_FACE_POS_X) {
172 u0 = uin.x;
173 } else if (user->identifiedInletBCFace == BC_FACE_NEG_Y ||
174 user->identifiedInletBCFace == BC_FACE_POS_Y) {
175 u0 = uin.y;
176 } else {
177 u0 = uin.z;
178 }
179 // now form the parabolic profile as before
180 normal_velocity_mag = 2.0 * u0 * (1.0 - 4.0 * r_sq);
181 }
182 }
183 break;
184 */
185 default:
186 LOG_ALLOW(LOCAL, LOG_WARNING, "Unrecognized FieldInitialization profile %d. Defaulting to zero.\n", simCtx->FieldInitialization);
187 normal_velocity_mag = 0.0;
188 break;
189 }
190
191 // Step B: Apply direction sign and set the correct contravariant component.
192 // The contravariant component U^n = v_n * Area_n, where v_n is the physical normal velocity.
193 if (normal_velocity_mag != 0.0) {
194 const PetscReal signed_normal_vel = normal_velocity_mag * flow_direction_sign*user->GridOrientation;
195
197 const PetscReal area_i = sqrt(csi_arr[k][j][i].x * csi_arr[k][j][i].x + csi_arr[k][j][i].y * csi_arr[k][j][i].y + csi_arr[k][j][i].z * csi_arr[k][j][i].z);
198
199 ucont_val.x = signed_normal_vel * area_i;
200
201 LOG_LOOP_ALLOW(GLOBAL,LOG_DEBUG,k,50," ucont_val.x = %.6f (signed_normal_vel=%.3f × area=%.4f)\n",ucont_val.x, signed_normal_vel, area_i);
202 }
204 const PetscReal area_j = sqrt(eta_arr[k][j][i].x * eta_arr[k][j][i].x + eta_arr[k][j][i].y * eta_arr[k][j][i].y + eta_arr[k][j][i].z * eta_arr[k][j][i].z);
205
206 ucont_val.y = signed_normal_vel * area_j;
207
208 LOG_LOOP_ALLOW(GLOBAL,LOG_DEBUG,k,50," ucont_val.y = %.6f (signed_normal_vel=%.3f × area=%.4f)\n",ucont_val.y, signed_normal_vel, area_j);
209 }
210 else { // Z-inlet
211 const PetscReal area_k = sqrt(zet_arr[k][j][i].x * zet_arr[k][j][i].x + zet_arr[k][j][i].y * zet_arr[k][j][i].y + zet_arr[k][j][i].z * zet_arr[k][j][i].z);
212
213 ucont_val.z = signed_normal_vel * area_k;
214
215 LOG_LOOP_ALLOW(GLOBAL,LOG_DEBUG,k,50," i,j,k,ucont_val.z = %d, %d, %d, %.6f (signed_normal_vel=%.3f × area=%.4f)\n",i,j,k,ucont_val.z, signed_normal_vel, area_k);
216 }
217 }
218 ucont_arr[k][j][i] = ucont_val;
219 } // end if(is_interior)
220 }
221 }
222 }
223 ierr = DMDAVecRestoreArray(fieldDM, fieldVec, &ucont_arr); CHKERRQ(ierr);
224
225 // --- 5. Cleanup: Restore arrays and free temporary memory ---
226 ierr = DMDAVecRestoreArrayRead(user->fda, localCoor, &coor_arr); CHKERRQ(ierr);
227 ierr = DMDAVecRestoreArrayRead(user->fda, user->lCsi, &csi_arr); CHKERRQ(ierr);
228 ierr = DMDAVecRestoreArrayRead(user->fda, user->lEta, &eta_arr); CHKERRQ(ierr);
229 ierr = DMDAVecRestoreArrayRead(user->fda, user->lZet, &zet_arr); CHKERRQ(ierr);
230
231 if (cent_coor) {
232 ierr = Deallocate3DArray(cent_coor, zm_cell, ym_cell); CHKERRQ(ierr);
233 }
234
235 PetscFunctionReturn(0);
236}
#define InterpolateFieldFromCornerToCenter(field, centfield, user)
Generic macro to call the appropriate interpolation function based on the field type.
#define LOG_LOOP_ALLOW(scope, level, iterVar, interval, fmt,...)
Logs a message inside a loop, but only every interval iterations.
Definition logging.h:319
#define LOCAL
Logging scope definitions for controlling message output.
Definition logging.h:44
#define GLOBAL
Scope for global logging across all processes.
Definition logging.h:45
#define LOG_ALLOW(scope, level, fmt,...)
Logging macro that checks both the log level and whether the calling function is in the allowed-funct...
Definition logging.h:207
@ LOG_INFO
Informational messages about program execution.
Definition logging.h:32
@ LOG_WARNING
Non-critical issues that warrant attention.
Definition logging.h:30
@ LOG_DEBUG
Detailed debugging information.
Definition logging.h:33
#define Allocate3DArray(array, nz, ny, nx)
Definition setup.h:28
#define Deallocate3DArray(array, nz, ny)
Definition setup.h:35
PetscErrorCode GetGhostedCellRange(const DMDALocalInfo *info_nodes, const RankNeighbors *neighbors, PetscInt dim, PetscInt *xs_cell_global_out, PetscInt *xm_cell_local_out)
Gets the global cell range for a rank, including boundary cells.
Definition setup.c:1452
BCFace identifiedInletBCFace
Definition variables.h:650
SimCtx * simCtx
Back-pointer to the master simulation context.
Definition variables.h:633
Vec lZet
Definition variables.h:673
PetscInt FieldInitialization
Definition variables.h:550
RankNeighbors neighbors
Definition variables.h:642
Vec Ucont
Definition variables.h:657
PetscScalar x
Definition variables.h:100
Vec lCsi
Definition variables.h:673
PetscScalar z
Definition variables.h:100
Cmpnts InitialConstantContra
Definition variables.h:551
PetscInt GridOrientation
Definition variables.h:643
PetscScalar y
Definition variables.h:100
Vec lEta
Definition variables.h:673
@ BC_FACE_NEG_X
Definition variables.h:200
@ BC_FACE_POS_Y
Definition variables.h:201
@ BC_FACE_POS_X
Definition variables.h:200
@ BC_FACE_NEG_Y
Definition variables.h:201
A 3D point or vector with PetscScalar components.
Definition variables.h:99
The master context for the entire simulation.
Definition variables.h:513
Here is the caller graph for this function:

◆ InitializeEulerianState()

PetscErrorCode InitializeEulerianState ( SimCtx simCtx)

High-level orchestrator to set the complete initial state of the Eulerian solver.

This function is called once from main() before the time loop begins. It inspects the simulation context to determine whether to perform a fresh start (t=0) or restart from saved files. It then delegates to the appropriate helper function. Finally, it initializes the solver's history vectors (Ucont_o, P_o, etc.) to ensure the first time step has the necessary data.

Definition at line 352 of file initialcondition.c.

353{
354 PetscErrorCode ierr;
355 UserCtx *user_finest = simCtx->usermg.mgctx[simCtx->usermg.mglevels - 1].user;
356
357 PetscFunctionBeginUser;
358
359 LOG_ALLOW(GLOBAL, LOG_INFO, "--- Initializing Eulerian State ---\n");
360
361 if (simCtx->rstart_flg) {
362 LOG_ALLOW(GLOBAL, LOG_INFO, "Starting from RESTART files (t=%.4f, step=%d).\n",
363 simCtx->StartTime, simCtx->StartStep);
364 ierr = SetInitialFluidState_Restart(simCtx); CHKERRQ(ierr);
365 } else {
366 LOG_ALLOW(GLOBAL, LOG_INFO, "Performing a FRESH START (t=0, step=0).\n");
367 ierr = SetInitialFluidState_FreshStart(simCtx); CHKERRQ(ierr);
368 }
369
370 // This crucial step, taken from the end of the legacy setup, ensures
371 // that the history vectors (Ucont_o, Ucont_rm1, etc.) are correctly
372 // populated before the first call to the time-stepping loop.
373 for (PetscInt bi = 0; bi < simCtx->block_number; bi++) {
374 ierr = UpdateSolverHistoryVectors(&user_finest[bi]); CHKERRQ(ierr);
375 }
376
377 LOG_ALLOW(GLOBAL, LOG_INFO, "--- Eulerian State Initialized and History Vectors Populated ---\n");
378 PetscFunctionReturn(0);
379}
static PetscErrorCode SetInitialFluidState_FreshStart(SimCtx *simCtx)
(HELPER) Sets the t=0 fluid state for all blocks.
static PetscErrorCode SetInitialFluidState_Restart(SimCtx *simCtx)
(HELPER) Reads fluid state for all blocks from restart files.
PetscErrorCode UpdateSolverHistoryVectors(UserCtx *user)
Copies the current time step's solution fields into history vectors (e.g., U(t_n) -> U_o,...
Definition simulation.c:23
UserCtx * user
Definition variables.h:418
PetscInt block_number
Definition variables.h:562
PetscReal StartTime
Definition variables.h:526
UserMG usermg
Definition variables.h:599
PetscBool rstart_flg
Definition variables.h:528
PetscInt StartStep
Definition variables.h:523
PetscInt mglevels
Definition variables.h:425
MGCtx * mgctx
Definition variables.h:428
User-defined context containing data specific to a single computational grid level.
Definition variables.h:630
Here is the call graph for this function:
Here is the caller graph for this function: