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

Implements the analytical solution engine for initializing or driving the simulation. More...

#include "AnalyticalSolutions.h"
#include <petscmath.h>
Include dependency graph for AnalyticalSolutions.c:

Go to the source code of this file.

Macros

#define __FUNCT__   "SetAnalyticalGridInfo"
 
#define __FUNCT__   "AnalyticalSolutionEngine"
 
#define __FUNCT__   "SetAnalyticalSolution_TGV3D"
 
#define __FUNCT__   "SetAnalyticalSolution_ZeroFlow"
 
#define __FUNCT__   "SetAnalyticalSolutionForParticles_TGV3D"
 
#define __FUNCT__   "SetAnalyticalSolutionForParticles"
 

Functions

static PetscErrorCode SetAnalyticalSolution_TGV3D (SimCtx *simCtx)
 Sets the non-dimensional velocity and pressure fields to the 3D Taylor-Green Vortex solution.
 
static PetscErrorCode SetAnalyticalSolution_ZeroFlow (SimCtx *simCtx)
 Sets all fluid fields to zero for a quiescent (zero-flow) background.
 
PetscBool AnalyticalTypeRequiresCustomGeometry (const char *analytical_type)
 Reports whether an analytical type requires custom geometry/decomposition logic.
 
PetscErrorCode SetAnalyticalGridInfo (UserCtx *user)
 Sets the grid domain and resolution for analytical solution cases.
 
PetscErrorCode AnalyticalSolutionEngine (SimCtx *simCtx)
 Dispatches to the appropriate analytical solution function based on simulation settings.
 
static PetscErrorCode SetAnalyticalSolutionForParticles_TGV3D (Vec tempVec, SimCtx *simCtx)
 Sets the TGV3D analytical velocity solution for particles.
 
PetscErrorCode SetAnalyticalSolutionForParticles (Vec tempVec, SimCtx *simCtx)
 Applies the analytical solution to particle velocity vector.
 

Detailed Description

Implements the analytical solution engine for initializing or driving the simulation.

This file provides a modular and extensible framework for applying analytical solutions to the Eulerian fields. The primary entry point is AnalyticalSolutionEngine, which acts as a dispatcher based on user configuration.

— DESIGN PHILOSOPHY —

  1. Non-Dimensional Core: All calculations within this engine are performed in non-dimensional units (e.g., reference velocity U_ref=1.0, reference length L_ref=1.0). This is critical for consistency with the core numerical solver, which also operates on non-dimensional equations. The simCtx->scaling parameters are intentionally NOT used here; they are reserved for dimensionalization during I/O and post-processing only.
  2. Separation of Concerns: The role of this engine is to set the physical state of the fluid at a given time t. This involves:
    • Setting the values of fields (Ucat, P) in the physical interior of the domain.
    • Declaring the physical values on the boundaries by populating the boundary condition vector (user->Bcs.Ubcs). It does NOT implement the numerical scheme for ghost cells. Instead, after setting the physical state, it relies on the solver's standard utility functions (UpdateDummyCells, UpdateCornerNodes) to correctly populate all ghost cell layers.
  3. Extensibility: The dispatcher design makes it straightforward to add new analytical solutions. A developer only needs to add a new else if condition and a corresponding static implementation function, without modifying any other part of the solver.

Definition in file AnalyticalSolutions.c.

Macro Definition Documentation

◆ __FUNCT__ [1/6]

#define __FUNCT__   "SetAnalyticalGridInfo"

Definition at line 45 of file AnalyticalSolutions.c.

◆ __FUNCT__ [2/6]

#define __FUNCT__   "AnalyticalSolutionEngine"

Definition at line 45 of file AnalyticalSolutions.c.

◆ __FUNCT__ [3/6]

#define __FUNCT__   "SetAnalyticalSolution_TGV3D"

Definition at line 45 of file AnalyticalSolutions.c.

◆ __FUNCT__ [4/6]

#define __FUNCT__   "SetAnalyticalSolution_ZeroFlow"

Definition at line 45 of file AnalyticalSolutions.c.

◆ __FUNCT__ [5/6]

#define __FUNCT__   "SetAnalyticalSolutionForParticles_TGV3D"

Definition at line 45 of file AnalyticalSolutions.c.

◆ __FUNCT__ [6/6]

#define __FUNCT__   "SetAnalyticalSolutionForParticles"

Definition at line 45 of file AnalyticalSolutions.c.

Function Documentation

◆ SetAnalyticalSolution_TGV3D()

static PetscErrorCode SetAnalyticalSolution_TGV3D ( SimCtx simCtx)
static

Sets the non-dimensional velocity and pressure fields to the 3D Taylor-Green Vortex solution.

This function implements the classical decaying Taylor-Green vortex in non-dimensional form, consistent with the solver's internal state. It is fully compatible with the solver's curvilinear, parallel, cell-centered architecture.

WORKFLOW
  1. Defines non-dimensional parameters for the TGV solution (V0=1.0, rho=1.0).
  2. Defines correct loop bounds (lxs, lxe, etc.) to iterate over ONLY the physical interior cells owned by the current MPI rank.
  3. Sets the interior cell-centered velocity (Ucat) and pressure (P) by evaluating the analytical formula at the cell-center coordinates (read from the local user->Cent vector).
  4. Sets the boundary condition vector (user->Bcs.Ubcs) by evaluating the analytical velocity at the true physical face-center locations (read from the local user->Centx, user->Centy, and user->Centz vectors).
  5. Manually sets the pressure ghost cells on the faces (excluding edges/corners) to enforce a zero-gradient (Neumann) boundary condition, which is standard for pressure.
  6. Calls the solver's utility functions UpdateDummyCells (for Ucat) and UpdateCornerNodes (for both Ucat and P) to finalize all ghost cell values, ensuring a fully consistent field.
Parameters
[in]simCtxThe main simulation context.
Returns
PetscErrorCode 0 on success, or a PETSc error code on failure.

Definition at line 250 of file AnalyticalSolutions.c.

251{
252 PetscErrorCode ierr;
253 UserCtx *user_finest = simCtx->usermg.mgctx[simCtx->usermg.mglevels - 1].user;
254
255 // --- NON-DIMENSIONAL TGV Parameters ---
256 const PetscReal V0 = 1.0; // Non-dimensional reference velocity.
257 const PetscReal rho = 1.0; // Non-dimensional reference density.
258 const PetscReal p0 = 0.0; // Non-dimensional reference pressure.
259
260 // Kinematic viscosity is derived from the non-dimensional Reynolds number.
261 const PetscReal nu = (simCtx->ren > 0) ? (1.0 / simCtx->ren) : 0.0;
262
263 const PetscReal k = 1.0; // Wavenumber, assumes a non-dimensional [0, 2*pi] domain.
264 const PetscReal t = simCtx->ti;
265
266 LOG_ALLOW(GLOBAL,LOG_TRACE,"TGV Setup: t = %.4f, V0* = %.4f, rho* = %.4f, k = %.4f, p0* = %4.f, nu = %.6f.\n",simCtx->ti,V0,rho,k,p0,nu);
267
268 const PetscReal vel_decay = exp(-2.0 * nu * k * k * t);
269 const PetscReal prs_decay = exp(-4.0 * nu * k * k * t);
270
271 PetscFunctionBeginUser;
272
273 for (PetscInt bi = 0; bi < simCtx->block_number; bi++) {
274 UserCtx* user = &user_finest[bi];
275 DMDALocalInfo info = user->info;
276 PetscInt xs = info.xs, xe = info.xs + info.xm;
277 PetscInt ys = info.ys, ye = info.ys + info.ym;
278 PetscInt zs = info.zs, ze = info.zs + info.zm;
279 PetscInt mx = info.mx, my = info.my, mz = info.mz;
280
281 Cmpnts ***ucat, ***ubcs;
282 const Cmpnts ***cent, ***cent_x, ***cent_y, ***cent_z;
283 PetscReal ***p;
284
285 // Define loop bounds for physical interior cells owned by this rank.
286 PetscInt lxs = (xs == 0) ? xs + 1 : xs, lxe = (xe == mx) ? xe - 1 : xe;
287 PetscInt lys = (ys == 0) ? ys + 1 : ys, lye = (ye == my) ? ye - 1 : ye;
288 PetscInt lzs = (zs == 0) ? zs + 1 : zs, lze = (ze == mz) ? ze - 1 : ze;
289
290 // --- Get Arrays ---
291 ierr = DMDAVecGetArray(user->fda, user->Ucat, &ucat); CHKERRQ(ierr);
292 ierr = DMDAVecGetArray(user->da, user->P, &p); CHKERRQ(ierr);
293 ierr = DMDAVecGetArray(user->fda, user->Bcs.Ubcs, &ubcs); CHKERRQ(ierr);
294 ierr = DMDAVecGetArrayRead(user->fda, user->Cent, &cent); CHKERRQ(ierr);
295 ierr = DMDAVecGetArrayRead(user->fda, user->Centx, &cent_x); CHKERRQ(ierr);
296 ierr = DMDAVecGetArrayRead(user->fda, user->Centy, &cent_y); CHKERRQ(ierr);
297 ierr = DMDAVecGetArrayRead(user->fda, user->Centz, &cent_z); CHKERRQ(ierr);
298
299 // --- Set INTERIOR cell-centered velocity (Ucat) ---
300 for (PetscInt k_cell = lzs; k_cell < lze; k_cell++) {
301 for (PetscInt j_cell = lys; j_cell < lye; j_cell++) {
302 for (PetscInt i_cell = lxs; i_cell < lxe; i_cell++) {
303 const PetscReal cx = cent[k_cell][j_cell][i_cell].x, cy = cent[k_cell][j_cell][i_cell].y, cz = cent[k_cell][j_cell][i_cell].z;
304 ucat[k_cell][j_cell][i_cell].x = V0 * sin(k*cx) * cos(k*cy) * cos(k*cz) * vel_decay;
305 ucat[k_cell][j_cell][i_cell].y = -V0 * cos(k*cx) * sin(k*cy) * cos(k*cz) * vel_decay;
306 ucat[k_cell][j_cell][i_cell].z = 0.0;
307 }
308 }
309 }
310
311
312 // --- Set INTERIOR cell-centered pressure (P) ---
313 for (PetscInt k_cell = lzs; k_cell < lze; k_cell++) {
314 for (PetscInt j_cell = lys; j_cell < lye; j_cell++) {
315 for (PetscInt i_cell = lxs; i_cell < lxe; i_cell++) {
316 const PetscReal cx = cent[k_cell][j_cell][i_cell].x, cy = cent[k_cell][j_cell][i_cell].y;
317 p[k_cell][j_cell][i_cell] = p0 + (rho * V0 * V0 / 4.0) * (cos(2*k*cx) + cos(2*k*cy)) * prs_decay;
318 }
319 }
320 }
321
322
323 // --- Set BOUNDARY condition vector for velocity (Ubcs) ---
324 if (xs == 0) for (PetscInt k=zs; k<ze; k++) for (PetscInt j=ys; j<ye; j++) {
325 const PetscReal fcx=cent_x[k][j][xs].x, fcy=cent_x[k][j][xs].y, fcz=cent_x[k][j][xs].z;
326 ubcs[k][j][xs].x = V0*sin(k*fcx)*cos(k*fcy)*cos(k*fcz)*vel_decay; ubcs[k][j][xs].y = -V0*cos(k*fcx)*sin(k*fcy)*cos(k*fcz)*vel_decay; ubcs[k][j][xs].z = 0.0;
327 }
328 if (xe == mx) for (PetscInt k=zs; k<ze; k++) for (PetscInt j=ys; j<ye; j++) {
329 const PetscReal fcx=cent_x[k][j][xe-1].x, fcy=cent_x[k][j][xe-1].y, fcz=cent_x[k][j][xe-1].z;
330 ubcs[k][j][xe-1].x = V0*sin(k*fcx)*cos(k*fcy)*cos(k*fcz)*vel_decay; ubcs[k][j][xe-1].y = -V0*cos(k*fcx)*sin(k*fcy)*cos(k*fcz)*vel_decay; ubcs[k][j][xe-1].z = 0.0;
331 }
332 if (ys == 0) for (PetscInt k=zs; k<ze; k++) for (PetscInt i=xs; i<xe; i++) {
333 const PetscReal fcx=cent_y[k][ys][i].x, fcy=cent_y[k][ys][i].y, fcz=cent_y[k][ys][i].z;
334 ubcs[k][ys][i].x = V0*sin(k*fcx)*cos(k*fcy)*cos(k*fcz)*vel_decay; ubcs[k][ys][i].y = -V0*cos(k*fcx)*sin(k*fcy)*cos(k*fcz)*vel_decay; ubcs[k][ys][i].z = 0.0;
335 }
336 if (ye == my) for (PetscInt k=zs; k<ze; k++) for (PetscInt i=xs; i<xe; i++) {
337 const PetscReal fcx=cent_y[k][ye-1][i].x, fcy=cent_y[k][ye-1][i].y, fcz=cent_y[k][ye-1][i].z;
338 ubcs[k][ye-1][i].x = V0*sin(k*fcx)*cos(k*fcy)*cos(k*fcz)*vel_decay; ubcs[k][ye-1][i].y = -V0*cos(k*fcx)*sin(k*fcy)*cos(k*fcz)*vel_decay; ubcs[k][ye-1][i].z = 0.0;
339 }
340 if (zs == 0) for (PetscInt j=ys; j<ye; j++) for (PetscInt i=xs; i<xe; i++) {
341 const PetscReal fcx=cent_z[zs][j][i].x, fcy=cent_z[zs][j][i].y, fcz=cent_z[zs][j][i].z;
342 ubcs[zs][j][i].x = V0*sin(k*fcx)*cos(k*fcy)*cos(k*fcz)*vel_decay; ubcs[zs][j][i].y = -V0*cos(k*fcx)*sin(k*fcy)*cos(k*fcz)*vel_decay; ubcs[zs][j][i].z = 0.0;
343 }
344 if (ze == mz) for (PetscInt j=ys; j<ye; j++) for (PetscInt i=xs; i<xe; i++) {
345 const PetscReal fcx=cent_z[ze-1][j][i].x, fcy=cent_z[ze-1][j][i].y, fcz=cent_z[ze-1][j][i].z;
346 ubcs[ze-1][j][i].x = V0*sin(k*fcx)*cos(k*fcy)*cos(k*fcz)*vel_decay; ubcs[ze-1][j][i].y = -V0*cos(k*fcx)*sin(k*fcy)*cos(k*fcz)*vel_decay; ubcs[ze-1][j][i].z = 0.0;
347 }
348
349 // --- Set PRESSURE GHOST CELLS (Neumann BC: P_ghost = P_interior) ---
350 if (xs == 0) for (PetscInt k=lzs; k<lze; k++) for (PetscInt j=lys; j<lye; j++) p[k][j][xs] = p[k][j][xs+1];
351 if (xe == mx) for (PetscInt k=lzs; k<lze; k++) for (PetscInt j=lys; j<lye; j++) p[k][j][xe-1] = p[k][j][xe-2];
352
353 if (ys == 0) for (PetscInt k=lzs; k<lze; k++) for (PetscInt i=lxs; i<lxe; i++) p[k][ys][i] = p[k][ys+1][i];
354 if (ye == my) for (PetscInt k=lzs; k<lze; k++) for (PetscInt i=lxs; i<lxe; i++) p[k][ye-1][i] = p[k][ye-2][i];
355
356 if (zs == 0) for (PetscInt j=lys; j<lye; j++) for (PetscInt i=lxs; i<lxe; i++) p[zs][j][i] = p[zs+1][j][i];
357 if (ze == mz) for (PetscInt j=lys; j<lye; j++) for (PetscInt i=lxs; i<lxe; i++) p[ze-1][j][i] = p[ze-2][j][i];
358
359 // --- Restore all arrays ---
360 ierr = DMDAVecRestoreArray(user->fda, user->Ucat, &ucat); CHKERRQ(ierr);
361 ierr = DMDAVecRestoreArray(user->da, user->P, &p); CHKERRQ(ierr);
362 ierr = DMDAVecRestoreArray(user->fda, user->Bcs.Ubcs, &ubcs); CHKERRQ(ierr);
363 ierr = DMDAVecRestoreArrayRead(user->fda, user->Cent, &cent); CHKERRQ(ierr);
364 ierr = DMDAVecRestoreArrayRead(user->fda, user->Centx, &cent_x); CHKERRQ(ierr);
365 ierr = DMDAVecRestoreArrayRead(user->fda, user->Centy, &cent_y); CHKERRQ(ierr);
366 ierr = DMDAVecRestoreArrayRead(user->fda, user->Centz, &cent_z); CHKERRQ(ierr);
367
368 // Pre-Dummy cell update synchronization.
369 ierr = UpdateLocalGhosts(user,"Ucat");
370 ierr = UpdateLocalGhosts(user,"P");
371
372 // --- Finalize all ghost cell values ---
373 ierr = UpdateDummyCells(user); CHKERRQ(ierr);
374 ierr = UpdateCornerNodes(user); CHKERRQ(ierr);
375
376 // Final Synchronization.
377 ierr = UpdateLocalGhosts(user,"Ucat");
378 ierr = UpdateLocalGhosts(user,"P");
379
380 }
381
382 PetscFunctionReturn(0);
383}
PetscErrorCode UpdateDummyCells(UserCtx *user)
Updates the dummy cells (ghost nodes) on the faces of the local domain for NON-PERIODIC boundaries.
PetscErrorCode UpdateCornerNodes(UserCtx *user)
Updates the corner and edge ghost nodes of the local domain by averaging.
#define GLOBAL
Scope for global logging across all processes.
Definition logging.h:46
#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:200
@ LOG_TRACE
Very fine-grained tracing information for in-depth debugging.
Definition logging.h:33
PetscErrorCode UpdateLocalGhosts(UserCtx *user, const char *fieldName)
Updates the local vector (including ghost points) from its corresponding global vector.
Definition setup.c:1162
UserCtx * user
Definition variables.h:474
PetscInt block_number
Definition variables.h:654
Vec Centz
Definition variables.h:782
UserMG usermg
Definition variables.h:703
PetscReal ren
Definition variables.h:636
Vec Ubcs
Boundary condition velocity values. (Comment: "An ugly hack, waste of memory")
Definition variables.h:121
PetscScalar x
Definition variables.h:101
Vec Centx
Definition variables.h:782
BCS Bcs
Definition variables.h:754
PetscScalar z
Definition variables.h:101
Vec Ucat
Definition variables.h:760
PetscInt mglevels
Definition variables.h:481
DMDALocalInfo info
Definition variables.h:740
PetscScalar y
Definition variables.h:101
Vec Cent
Definition variables.h:781
MGCtx * mgctx
Definition variables.h:484
Vec Centy
Definition variables.h:782
PetscReal ti
Definition variables.h:598
A 3D point or vector with PetscScalar components.
Definition variables.h:100
User-defined context containing data specific to a single computational grid level.
Definition variables.h:733
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetAnalyticalSolution_ZeroFlow()

static PetscErrorCode SetAnalyticalSolution_ZeroFlow ( SimCtx simCtx)
static

Sets all fluid fields to zero for a quiescent (zero-flow) background.

Used for Brownian-motion validation: particles diffuse from a point source with no advection, so the carrier flow is identically zero everywhere. The ghost-cell sequence is identical to that used by TGV3D, ensuring a fully consistent field state before the first time step.

Parameters
[in]simCtxThe main simulation context.
Returns
PetscErrorCode 0 on success.

Definition at line 398 of file AnalyticalSolutions.c.

399{
400 PetscErrorCode ierr;
401 UserCtx *user_finest = simCtx->usermg.mgctx[simCtx->usermg.mglevels - 1].user;
402
403 PetscFunctionBeginUser;
404
405 for (PetscInt bi = 0; bi < simCtx->block_number; bi++) {
406 UserCtx *user = &user_finest[bi];
407
408 ierr = VecZeroEntries(user->Ucat); CHKERRQ(ierr);
409 ierr = VecZeroEntries(user->P); CHKERRQ(ierr);
410 ierr = VecZeroEntries(user->Bcs.Ubcs); CHKERRQ(ierr);
411
412 // Ghost-cell finalization — identical sequence to TGV3D
413 ierr = UpdateLocalGhosts(user, "Ucat"); CHKERRQ(ierr);
414 ierr = UpdateLocalGhosts(user, "P"); CHKERRQ(ierr);
415 ierr = UpdateDummyCells(user); CHKERRQ(ierr);
416 ierr = UpdateCornerNodes(user); CHKERRQ(ierr);
417 ierr = UpdateLocalGhosts(user, "Ucat"); CHKERRQ(ierr);
418 ierr = UpdateLocalGhosts(user, "P"); CHKERRQ(ierr);
419 }
420
421 PetscFunctionReturn(0);
422}
Here is the call graph for this function:
Here is the caller graph for this function:

◆ AnalyticalTypeRequiresCustomGeometry()

PetscBool AnalyticalTypeRequiresCustomGeometry ( const char *  analytical_type)

Reports whether an analytical type requires custom geometry/decomposition logic.

Analytical types returning PETSC_TRUE are expected to route through SetAnalyticalGridInfo. Types returning PETSC_FALSE should use the standard programmatic grid parser fallback.

Parameters
analytical_typeAnalytical solution type string.
Returns
PETSC_TRUE if custom geometry is required, PETSC_FALSE otherwise.

Definition at line 38 of file AnalyticalSolutions.c.

39{
40 if (!analytical_type) return PETSC_FALSE;
41 return (strcmp(analytical_type, "TGV3D") == 0) ? PETSC_TRUE : PETSC_FALSE;
42}
Here is the caller graph for this function:

◆ SetAnalyticalGridInfo()

PetscErrorCode SetAnalyticalGridInfo ( UserCtx user)

Sets the grid domain and resolution for analytical solution cases.

This function is called when eulerianSource is "analytical". It is responsible for automatically configuring the grid based on the chosen AnalyticalSolutionType.

TGV3D Multi-Block Decomposition
If the analytical solution is "TGV3D", this function automatically decomposes the required [0, 2*PI] physical domain among the available blocks.
  • **Single Block (nblk=1):** The single block is assigned the full [0, 2*PI] domain.
  • **Multiple Blocks (nblk>1):** It requires that the number of blocks be a perfect square (e.g., 4, 9, 16). It then arranges the blocks in a sqrt(nblk) by sqrt(nblk) grid in the X-Y plane, partitioning the [0, 2*PI] domain in X and Y accordingly. The Z domain for all blocks remains [0, 2*PI]. If nblk is not a perfect square, the simulation is aborted with an error.

Grid resolution (IM/JM/KM) is expected to be pre-populated in user before this function is called.

Parameters
userPointer to the UserCtx for a specific block. The function will populate the geometric fields (IM, JM, KM, Min_X, Max_X, etc.) within this struct.
Returns
PetscErrorCode 0 on success, or a PETSc error code on failure.

Definition at line 70 of file AnalyticalSolutions.c.

71{
72 SimCtx *simCtx = user->simCtx;
73 PetscInt nblk = simCtx->block_number;
74 PetscInt block_index = user->_this;
75
76 PetscFunctionBeginUser;
78
80 SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONGSTATE,
81 "SetAnalyticalGridInfo called for analytical type '%s' that does not require custom geometry.",
83 }
84 if (user->IM <= 0 || user->JM <= 0 || user->KM <= 0) {
85 SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONGSTATE,
86 "Analytical grid resolution is not initialized. Ensure IM/JM/KM are preloaded before SetAnalyticalGridInfo.");
87 }
88
89 if (strcmp(simCtx->AnalyticalSolutionType, "TGV3D") == 0) {
90 LOG_ALLOW_SYNC(GLOBAL, LOG_DEBUG, "Rank %d: Configuring grid for TGV3D analytical solution, block %d.\n", simCtx->rank, block_index);
91
92 if (nblk == 1) {
93 // --- Single Block Case ---
94 if (block_index == 0) {
95 LOG_ALLOW(GLOBAL, LOG_INFO, "Single block detected. Setting domain to [0, 2*PI].\n");
96 }
97 user->Min_X = 0.0; user->Max_X = 2.0 * PETSC_PI;
98 user->Min_Y = 0.0; user->Max_Y = 2.0 * PETSC_PI;
99 user->Min_Z = 0.0; user->Max_Z = 0.2 * PETSC_PI; //2.0 * PETSC_PI;
100
101 } else { // --- Multi-Block Case ---
102 PetscReal s = sqrt((PetscReal)nblk);
103
104 // Validate that nblk is a perfect square.
105 if (fabs(s - floor(s)) > 1e-9) {
106 SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_INCOMP,
107 "\n\n*** CONFIGURATION ERROR FOR TGV3D ***\n"
108 "For multi-block TGV3D cases, the number of blocks must be a perfect square (e.g., 4, 9, 16).\n"
109 "You have specified %d blocks. Please adjust `-block_number`.\n", nblk);
110 }
111 PetscInt blocks_per_dim = (PetscInt)s;
112
113 if (block_index == 0) {
114 LOG_ALLOW(GLOBAL, LOG_INFO, "%d blocks detected. Decomposing domain into a %d x %d grid in the X-Y plane.\n", nblk, blocks_per_dim, blocks_per_dim);
115 }
116
117 // Determine the (row, col) position of this block in the 2D decomposition
118 PetscInt row = block_index / blocks_per_dim;
119 PetscInt col = block_index % blocks_per_dim;
120
121 // Calculate the width/height of each sub-domain
122 PetscReal block_width = (2.0 * PETSC_PI) / (PetscReal)blocks_per_dim;
123 PetscReal block_height = (2.0 * PETSC_PI) / (PetscReal)blocks_per_dim;
124
125 // Assign this block its specific sub-domain
126 user->Min_X = col * block_width;
127 user->Max_X = (col + 1) * block_width;
128 user->Min_Y = row * block_height;
129 user->Max_Y = (row + 1) * block_height;
130 user->Min_Z = 0.0;
131 user->Max_Z = 2.0 * PETSC_PI; // Z-domain is not decomposed
132 }
133 }
134 /*
135 * --- EXTENSIBILITY HOOK ---
136 * To add another analytical case with special grid requirements:
137 *
138 * else if (strcmp(simCtx->AnalyticalSolutionType, "ChannelFlow") == 0) {
139 * // ... implement logic to set domain for ChannelFlow case ...
140 * }
141 */
142 else {
143 SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_UNKNOWN_TYPE,
144 "Analytical type '%s' has no custom geometry implementation.",
145 simCtx->AnalyticalSolutionType);
146 }
147
148 // We can also read stretching ratios, as they are independent of the domain size
149 // For simplicity, we assume uniform grid unless specified.
150 user->rx = 1.0; user->ry = 1.0; user->rz = 1.0;
151
152 LOG_ALLOW(LOCAL, LOG_DEBUG, "Rank %d: Block %d grid resolution set: IM=%d, JM=%d, KM=%d\n",
153 simCtx->rank, block_index, user->IM, user->JM, user->KM);
154 LOG_ALLOW(LOCAL, LOG_DEBUG, "Rank %d: Block %d final bounds: X=[%.4f, %.4f], Y=[%.4f, %.4f], Z=[%.4f, %.4f]\n",
155 simCtx->rank, block_index, user->Min_X, user->Max_X, user->Min_Y, user->Max_Y, user->Min_Z, user->Max_Z);
156
158 PetscFunctionReturn(0);
159}
PetscBool AnalyticalTypeRequiresCustomGeometry(const char *analytical_type)
Reports whether an analytical type requires custom geometry/decomposition logic.
#define LOG_ALLOW_SYNC(scope, level, fmt,...)
----— DEBUG ---------------------------------------— #define LOG_ALLOW(scope, level,...
Definition logging.h:267
#define LOCAL
Logging scope definitions for controlling message output.
Definition logging.h:45
#define PROFILE_FUNCTION_END
Marks the end of a profiled code block.
Definition logging.h:746
@ LOG_INFO
Informational messages about program execution.
Definition logging.h:31
@ LOG_DEBUG
Detailed debugging information.
Definition logging.h:32
#define PROFILE_FUNCTION_BEGIN
Marks the beginning of a profiled code block (typically a function).
Definition logging.h:737
PetscMPIInt rank
Definition variables.h:592
SimCtx * simCtx
Back-pointer to the master simulation context.
Definition variables.h:736
PetscReal Min_X
Definition variables.h:743
PetscInt KM
Definition variables.h:742
PetscInt _this
Definition variables.h:746
PetscReal ry
Definition variables.h:747
PetscReal Max_Y
Definition variables.h:743
PetscReal rz
Definition variables.h:747
PetscInt JM
Definition variables.h:742
PetscReal Min_Z
Definition variables.h:743
char AnalyticalSolutionType[PETSC_MAX_PATH_LEN]
Definition variables.h:621
PetscReal Max_X
Definition variables.h:743
PetscReal Min_Y
Definition variables.h:743
PetscInt IM
Definition variables.h:742
PetscReal rx
Definition variables.h:747
PetscReal Max_Z
Definition variables.h:743
The master context for the entire simulation.
Definition variables.h:589
Here is the call graph for this function:
Here is the caller graph for this function:

◆ AnalyticalSolutionEngine()

PetscErrorCode AnalyticalSolutionEngine ( SimCtx simCtx)

Dispatches to the appropriate analytical solution function based on simulation settings.

This function acts as a router. It reads the AnalyticalSolutionType string from the simulation context and calls the corresponding private implementation function (e.g., for Taylor-Green Vortex). This design keeps the main simulation code clean and makes it easy to add new analytical test cases.

Parameters
[in]simCtxThe main simulation context, containing all configuration and state.
Returns
PetscErrorCode 0 on success, or a PETSc error code on failure.

Definition at line 179 of file AnalyticalSolutions.c.

180{
181 PetscErrorCode ierr;
182 PetscFunctionBeginUser;
184
185 // -- Before any operation, here is defensive test to ensure that the Corner->Center Interpolation method works
186 //ierr = TestCornerToCenterInterpolation(&(simCtx->usermg.mgctx[simCtx->usermg.mglevels - 1]->user[0]));
187
188 // --- Dispatch based on the string provided by the user ---
189 if (strcmp(simCtx->AnalyticalSolutionType, "TGV3D") == 0) {
190 LOG_ALLOW(GLOBAL, LOG_DEBUG, "Applying Analytical Solution: 3D Taylor-Green Vortex (TGV3D).\n");
191 ierr = SetAnalyticalSolution_TGV3D(simCtx); CHKERRQ(ierr);
192 }
193 else if (strcmp(simCtx->AnalyticalSolutionType, "ZERO_FLOW") == 0) {
194 LOG_ALLOW(GLOBAL, LOG_DEBUG, "Applying Analytical Solution: Zero Background Flow (ZERO_FLOW).\n");
195 ierr = SetAnalyticalSolution_ZeroFlow(simCtx); CHKERRQ(ierr);
196 }
197 /*
198 * --- EXTENSIBILITY HOOK ---
199 * To add a new analytical solution (e.g., "ChannelFlow"):
200 * 1. Add an `else if` block here:
201 *
202 * else if (strcmp(simCtx->AnalyticalSolutionType, "ChannelFlow") == 0) {
203 * LOG_ALLOW(GLOBAL, LOG_DEBUG, "Applying Analytical Solution: Channel Flow.\n");
204 * ierr = SetAnalyticalSolution_ChannelFlow(simCtx); CHKERRQ(ierr);
205 * }
206 *
207 * 2. Implement the static function `SetAnalyticalSolution_ChannelFlow(SimCtx *simCtx)`
208 * below, following the TGV3D pattern.
209 */
210 else {
211 // If the type is unknown, raise a fatal error to prevent silent failures.
212 SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown AnalyticalSolutionType specified: '%s'", simCtx->AnalyticalSolutionType);
213 }
214
216 PetscFunctionReturn(0);
217}
static PetscErrorCode SetAnalyticalSolution_ZeroFlow(SimCtx *simCtx)
Sets all fluid fields to zero for a quiescent (zero-flow) background.
static PetscErrorCode SetAnalyticalSolution_TGV3D(SimCtx *simCtx)
Sets the non-dimensional velocity and pressure fields to the 3D Taylor-Green Vortex solution.
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetAnalyticalSolutionForParticles_TGV3D()

static PetscErrorCode SetAnalyticalSolutionForParticles_TGV3D ( Vec  tempVec,
SimCtx simCtx 
)
static

Sets the TGV3D analytical velocity solution for particles.

Computes the 3D Taylor-Green Vortex velocity field at each particle position. Assumes the vector contains interleaved xyz components [x0,y0,z0, x1,y1,z1, ...].

Parameters
tempVecThe PETSc Vec containing particle positions, will be overwritten with velocities.
simCtxThe simulation context containing time and Reynolds number.
Returns
PetscErrorCode Returns 0 on success.

Definition at line 436 of file AnalyticalSolutions.c.

437{
438 PetscErrorCode ierr;
439 PetscInt nLocal;
440 PetscReal *data;
441
442 PetscFunctionBeginUser;
443
444 // TGV3D parameters (matching your Eulerian implementation)
445 const PetscReal V0 = 1.0;
446 const PetscReal k = 1.0;
447 const PetscReal nu = (simCtx->ren > 0) ? (1.0 / simCtx->ren) : 0.0;
448 const PetscReal t = simCtx->ti;
449 const PetscReal vel_decay = exp(-2.0 * nu * k * k * t);
450
451 LOG_ALLOW(GLOBAL, LOG_DEBUG, "TGV3D Particles: t=%.4f, V0=%.4f, k=%.4f, nu=%.6f\n", t, V0, k, nu);
452
453 ierr = VecGetLocalSize(tempVec, &nLocal); CHKERRQ(ierr);
454 ierr = VecGetArray(tempVec, &data); CHKERRQ(ierr);
455
456 // Process particles: data is interleaved [x0,y0,z0, x1,y1,z1, ...]
457 for (PetscInt i = 0; i < nLocal; i += 3) {
458 const PetscReal x = data[i];
459 const PetscReal y = data[i+1];
460 const PetscReal z = data[i+2];
461
462 // TGV3D velocity field
463 data[i] = V0 * sin(k*x) * cos(k*y) * cos(k*z) * vel_decay; // u
464 data[i+1] = -V0 * cos(k*x) * sin(k*y) * cos(k*z) * vel_decay; // v
465 data[i+2] = 0.0; // w
466 }
467
468 ierr = VecRestoreArray(tempVec, &data); CHKERRQ(ierr);
469
470 PetscFunctionReturn(0);
471}
Here is the caller graph for this function:

◆ SetAnalyticalSolutionForParticles()

PetscErrorCode SetAnalyticalSolutionForParticles ( Vec  tempVec,
SimCtx simCtx 
)

Applies the analytical solution to particle velocity vector.

Dispatcher function that calls the appropriate analytical solution based on simCtx->AnalyticalSolutionType. Supports multiple solution types.

Parameters
tempVecThe PETSc Vec containing particle positions which will be used to store velocities.
simCtxThe simulation context.
Returns
PetscErrorCode Returns 0 on success.

Definition at line 485 of file AnalyticalSolutions.c.

486{
487 PetscErrorCode ierr;
488 PetscInt nLocal;
489 PetscReal *vels;
490
491 PetscFunctionBeginUser;
492
493 LOG_ALLOW(GLOBAL, LOG_DEBUG, "Type: %s\n",
494 simCtx->AnalyticalSolutionType ? simCtx->AnalyticalSolutionType : "default");
495
496 // Check for specific analytical solution types
497 if (simCtx->AnalyticalSolutionType && strcmp(simCtx->AnalyticalSolutionType, "TGV3D") == 0) {
498 LOG_ALLOW(GLOBAL, LOG_DEBUG, "Using TGV3D solution.\n");
499 ierr = SetAnalyticalSolutionForParticles_TGV3D(tempVec, simCtx); CHKERRQ(ierr);
500 return 0;
501 }
502
503 ierr = VecRestoreArray(tempVec, &vels); CHKERRQ(ierr);
504
505 PetscFunctionReturn(0);
506}
static PetscErrorCode SetAnalyticalSolutionForParticles_TGV3D(Vec tempVec, SimCtx *simCtx)
Sets the TGV3D analytical velocity solution for particles.
Here is the call graph for this function:
Here is the caller graph for this function: