PICurv 0.1.0
A Parallel Particle-In-Cell Solver for Curvilinear LES
Loading...
Searching...
No Matches
initialcondition.c File Reference

// Setup the Initial conditions for different cases. More...

#include "initialcondition.h"
Include dependency graph for initialcondition.c:

Go to the source code of this file.

Macros

#define __FUNCT__   "SetInitialInteriorField"
 
#define __FUNCT__   "FinalizeBlockState"
 
#define __FUNCT__   "SetInitialFluidState_FreshStart"
 
#define __FUNCT__   "SetInitialFluidState_Load"
 
#define __FUNCT__   "InitializeEulerianState"
 

Functions

PetscErrorCode SetInitialInteriorField (UserCtx *user, const char *fieldName)
 Sets the initial values for the INTERIOR of a specified Eulerian field.
 
static PetscErrorCode FinalizeBlockState (UserCtx *user)
 (HELPER) Finalizes a block's state by applying BCs, converting to Cartesian, and synchronizing ghost regions.
 
static PetscErrorCode SetInitialFluidState_FreshStart (SimCtx *simCtx)
 (HELPER) Sets the t=0 fluid state for all blocks.
 
static PetscErrorCode SetInitialFluidState_Load (SimCtx *simCtx)
 (HELPER) Reads fluid state for all blocks from restart files.
 
PetscErrorCode InitializeEulerianState (SimCtx *simCtx)
 High-level orchestrator to set the complete initial state of the Eulerian solver.
 

Detailed Description

// Setup the Initial conditions for different cases.

Test program for DMSwarm interpolation using the fdf-curvIB method.

Definition in file initialcondition.c.

Macro Definition Documentation

◆ __FUNCT__ [1/5]

#define __FUNCT__   "SetInitialInteriorField"

Definition at line 9 of file initialcondition.c.

◆ __FUNCT__ [2/5]

#define __FUNCT__   "FinalizeBlockState"

Definition at line 9 of file initialcondition.c.

◆ __FUNCT__ [3/5]

#define __FUNCT__   "SetInitialFluidState_FreshStart"

Definition at line 9 of file initialcondition.c.

◆ __FUNCT__ [4/5]

#define __FUNCT__   "SetInitialFluidState_Load"

Definition at line 9 of file initialcondition.c.

◆ __FUNCT__ [5/5]

#define __FUNCT__   "InitializeEulerianState"

Definition at line 9 of file initialcondition.c.

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 32 of file initialcondition.c.

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

◆ FinalizeBlockState()

static PetscErrorCode FinalizeBlockState ( UserCtx user)
static

(HELPER) Finalizes a block's state by applying BCs, converting to Cartesian, and synchronizing ghost regions.

Definition at line 253 of file initialcondition.c.

254{
255 PetscErrorCode ierr;
256 PetscFunctionBeginUser;
257
259
260 // This sequence ensures a fully consistent state for a single block.
261 // 1. Apply BCs using the information from InflowFlux/OutflowFlux.
262 ierr = FormBCS(user); CHKERRQ(ierr);
263 LOG_ALLOW(GLOBAL,LOG_TRACE," Boundary condition applied.\n");
264 // 2. Sync contravariant velocity field.
265 ierr = UpdateLocalGhosts(user, "Ucont"); CHKERRQ(ierr);
266 LOG_ALLOW(GLOBAL,LOG_TRACE," Ucont field ghosts updated.\n");
267
268 // 3. Convert to Cartesian velocity.
269 ierr = Contra2Cart(user); CHKERRQ(ierr);
270 LOG_ALLOW(GLOBAL,LOG_TRACE," Converted Ucont to Ucat.\n");
271
272 // 4. Sync the new Cartesian velocity field.
273 ierr = UpdateLocalGhosts(user, "Ucat"); CHKERRQ(ierr);
274 LOG_ALLOW(GLOBAL,LOG_TRACE," Ucat field ghosts updated.\n");
275
277
278 PetscFunctionReturn(0);
279}
PetscErrorCode FormBCS(UserCtx *user)
@ LOG_TRACE
Very fine-grained tracing information for in-depth debugging.
Definition logging.h:34
PetscErrorCode Contra2Cart(UserCtx *user)
Reconstructs Cartesian velocity (Ucat) at cell centers from contravariant velocity (Ucont) defined on...
Definition setup.c:2084
PetscErrorCode UpdateLocalGhosts(UserCtx *user, const char *fieldName)
Updates the local vector (including ghost points) from its corresponding global vector.
Definition setup.c:1091
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetInitialFluidState_FreshStart()

static PetscErrorCode SetInitialFluidState_FreshStart ( SimCtx simCtx)
static

(HELPER) Sets the t=0 fluid state for all blocks.

This function replicates the initialization sequence for a fresh start from the legacy code. It sets an initial guess for the interior, establishes the boundary conditions, and ensures the final state is consistent across MPI ranks.

Definition at line 291 of file initialcondition.c.

292{
293 PetscErrorCode ierr;
294 UserCtx *user_finest = simCtx->usermg.mgctx[simCtx->usermg.mglevels - 1].user;
295 PetscFunctionBeginUser;
296
298
299 for (PetscInt bi = 0; bi < simCtx->block_number; bi++) {
300 LOG_ALLOW(LOCAL, LOG_INFO, "Rank %d, Block %d: Setting t=0 state.\n", simCtx->rank, bi);
301
302 // 1. Set an initial guess for the INTERIOR of the domain.
303 // Replaces the legacy `if(InitialGuessOne)` block.
304
305 LOG_ALLOW(GLOBAL,LOG_TRACE," Initializing Interior Ucont field.\n");
306 ierr = SetInitialInteriorField(&user_finest[bi], "Ucont"); CHKERRQ(ierr);
307 LOG_ALLOW(GLOBAL,LOG_TRACE," Interior Ucont field initialized.\n");
308
309 // 2. Establish the initial boundary flux values and set inlet profiles.
310 // This sequence is critical and comes directly from the legacy setup.
311 LOG_ALLOW(GLOBAL,LOG_TRACE," Boundary condition Preparation steps initiated.\n");
312 ierr = InflowFlux(&user_finest[bi]); CHKERRQ(ierr);
313 ierr = OutflowFlux(&user_finest[bi]); CHKERRQ(ierr);
314 LOG_ALLOW(GLOBAL,LOG_TRACE," Boundary condition Preparation steps completed.\n");
315 // 3. Apply all boundary conditions, convert to Cartesian, and sync ghosts.
316 LOG_ALLOW(GLOBAL,LOG_TRACE," Boundary condition application and state finalization initiated.\n");
317 ierr = FinalizeBlockState(&user_finest[bi]); CHKERRQ(ierr);
318 LOG_ALLOW(GLOBAL,LOG_TRACE," Boundary condition application and state finalization complete.\n");
319 }
320
321 // If using multiple grid blocks, handle the interface conditions between them.
322 if (simCtx->block_number > 1) {
323 // LOG_ALLOW(GLOBAL, LOG_INFO, "Updating multi-block interfaces for t=0.\n");
324 // ierr = Block_Interface_U(user_finest); CHKERRQ(ierr);
325 // After interface update, ghost regions might be stale. Refresh them.
326 for (PetscInt bi = 0; bi < simCtx->block_number; bi++) {
327 ierr = UpdateLocalGhosts(&user_finest[bi], "Ucont"); CHKERRQ(ierr);
328 ierr = UpdateLocalGhosts(&user_finest[bi], "Ucat"); CHKERRQ(ierr);
329 }
330 }
331
333
334 PetscFunctionReturn(0);
335}
PetscErrorCode OutflowFlux(UserCtx *user)
Calculates the total outgoing flux through all OUTLET faces for reporting.
Definition Boundaries.c:932
PetscErrorCode InflowFlux(UserCtx *user)
Applies inlet boundary conditions based on the modern BC handling system.
Definition Boundaries.c:706
PetscErrorCode SetInitialInteriorField(UserCtx *user, const char *fieldName)
Sets the initial values for the INTERIOR of a specified Eulerian field.
static PetscErrorCode FinalizeBlockState(UserCtx *user)
(HELPER) Finalizes a block's state by applying BCs, converting to Cartesian, and synchronizing ghost ...
UserCtx * user
Definition variables.h:427
PetscMPIInt rank
Definition variables.h:541
PetscInt block_number
Definition variables.h:593
UserMG usermg
Definition variables.h:631
PetscInt mglevels
Definition variables.h:434
MGCtx * mgctx
Definition variables.h:437
User-defined context containing data specific to a single computational grid level.
Definition variables.h:661
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetInitialFluidState_Load()

static PetscErrorCode SetInitialFluidState_Load ( SimCtx simCtx)
static

(HELPER) Reads fluid state for all blocks from restart files.

Definition at line 342 of file initialcondition.c.

343{
344 PetscErrorCode ierr;
345 UserCtx *user_finest = simCtx->usermg.mgctx[simCtx->mglevels - 1].user;
346 PetscFunctionBeginUser;
347
349
350 for (PetscInt bi = 0; bi < simCtx->block_number; bi++) {
351 LOG_ALLOW(LOCAL, LOG_INFO, "Rank %d, Block %d: Reading restart files for step %d.\n",
352 simCtx->rank, bi, simCtx->StartStep);
353
354 // ReadSimulationFields handles all file I/O for one block.
355 ierr = ReadSimulationFields(&user_finest[bi], simCtx->StartStep); CHKERRQ(ierr);
356
357 // After reading from a file, the local ghost regions MUST be updated
358 // to ensure consistency across process boundaries for the first time step.
359 ierr = UpdateLocalGhosts(&user_finest[bi], "Ucont"); CHKERRQ(ierr);
360 ierr = UpdateLocalGhosts(&user_finest[bi], "Ucat"); CHKERRQ(ierr);
361 ierr = UpdateLocalGhosts(&user_finest[bi], "P"); CHKERRQ(ierr);
362 // ... add ghost updates for any other fields read from file ...
363 }
364
366
367 PetscFunctionReturn(0);
368}
PetscErrorCode ReadSimulationFields(UserCtx *user, PetscInt ti)
Reads binary field data for velocity, pressure, and other required vectors.
Definition io.c:1024
PetscInt StartStep
Definition variables.h:548
PetscInt mglevels
Definition variables.h:577
Here is the call graph for this function:
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 381 of file initialcondition.c.

382{
383 PetscErrorCode ierr;
384 UserCtx *user_finest = simCtx->usermg.mgctx[simCtx->usermg.mglevels - 1].user;
385
386 PetscFunctionBeginUser;
387
389
390 LOG_ALLOW(GLOBAL, LOG_INFO, "--- Initializing Eulerian State ---\n");
391
392 if (simCtx->StartStep > 0) {
393 LOG_ALLOW(GLOBAL, LOG_INFO, "Starting from RESTART files (t=%.4f, step=%d).\n",
394 simCtx->StartTime, simCtx->StartStep);
395 ierr = SetInitialFluidState_Load(simCtx); CHKERRQ(ierr);
396 } else {
397 LOG_ALLOW(GLOBAL, LOG_INFO, "Performing a FRESH START (t=0, step=0).\n");
398 if(strcmp(simCtx->eulerianSource,"solve")==0){
399 ierr = SetInitialFluidState_FreshStart(simCtx); CHKERRQ(ierr);
400 }else if(strcmp(simCtx->eulerianSource,"load")==0){
401
402 LOG_ALLOW(GLOBAL,LOG_INFO,"FRESH START in LOAD mode. Reading files (t=%.4f,step=%d).\n",
403 simCtx->StartTime,simCtx->StartStep);
404
405 ierr=SetInitialFluidState_Load(simCtx);CHKERRQ(ierr);
406 }
407 }
408
409 // This crucial step, taken from the end of the legacy setup, ensures
410 // that the history vectors (Ucont_o, Ucont_rm1, etc.) are correctly
411 // populated before the first call to the time-stepping loop.
412 for (PetscInt bi = 0; bi < simCtx->block_number; bi++) {
413 ierr = UpdateSolverHistoryVectors(&user_finest[bi]); CHKERRQ(ierr);
414 }
415
416 LOG_ALLOW(GLOBAL, LOG_INFO, "--- Eulerian State Initialized and History Vectors Populated ---\n");
417
419
420 PetscFunctionReturn(0);
421}
static PetscErrorCode SetInitialFluidState_Load(SimCtx *simCtx)
(HELPER) Reads fluid state for all blocks from restart files.
static PetscErrorCode SetInitialFluidState_FreshStart(SimCtx *simCtx)
(HELPER) Sets the t=0 fluid state for all blocks.
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
char eulerianSource[64]
Definition variables.h:557
PetscReal StartTime
Definition variables.h:551
Here is the call graph for this function:
Here is the caller graph for this function: