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

Go to the source code of this file.

Macros

#define __FUNCT__   "AnalyticalSolutionEngine"
 

Functions

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.
 
PetscErrorCode SetAnalyticalSolutionForParticles (Vec tempVec, SimCtx *simCtx)
 Applies the analytical solution to particle velocity vector.
 

Macro Definition Documentation

◆ __FUNCT__

#define __FUNCT__   "AnalyticalSolutionEngine"

Definition at line 49 of file AnalyticalSolutions.h.

Function Documentation

◆ 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.

After setting the domain bounds, it proceeds to read the grid resolution options (-im, -jm, -km) from the command line for the specific block.

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 63 of file AnalyticalSolutions.c.

64{
65 PetscErrorCode ierr;
66 SimCtx *simCtx = user->simCtx;
67 PetscInt nblk = simCtx->block_number;
68 PetscInt block_index = user->_this;
69
70 PetscFunctionBeginUser;
72
73 if (strcmp(simCtx->AnalyticalSolutionType, "TGV3D") == 0) {
74 LOG_ALLOW_SYNC(GLOBAL, LOG_DEBUG, "Rank %d: Configuring grid for TGV3D analytical solution, block %d.\n", simCtx->rank, block_index);
75
76 if (nblk == 1) {
77 // --- Single Block Case ---
78 if (block_index == 0) {
79 LOG_ALLOW(GLOBAL, LOG_INFO, "Single block detected. Setting domain to [0, 2*PI].\n");
80 }
81 user->Min_X = 0.0; user->Max_X = 2.0 * PETSC_PI;
82 user->Min_Y = 0.0; user->Max_Y = 2.0 * PETSC_PI;
83 user->Min_Z = 0.0; user->Max_Z = 0.2 * PETSC_PI; //2.0 * PETSC_PI;
84
85 } else { // --- Multi-Block Case ---
86 PetscReal s = sqrt((PetscReal)nblk);
87
88 // Validate that nblk is a perfect square.
89 if (fabs(s - floor(s)) > 1e-9) {
90 SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_INCOMP,
91 "\n\n*** CONFIGURATION ERROR FOR TGV3D ***\n"
92 "For multi-block TGV3D cases, the number of blocks must be a perfect square (e.g., 4, 9, 16).\n"
93 "You have specified %d blocks. Please adjust `-block_number`.\n", nblk);
94 }
95 PetscInt blocks_per_dim = (PetscInt)s;
96
97 if (block_index == 0) {
98 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);
99 }
100
101 // Determine the (row, col) position of this block in the 2D decomposition
102 PetscInt row = block_index / blocks_per_dim;
103 PetscInt col = block_index % blocks_per_dim;
104
105 // Calculate the width/height of each sub-domain
106 PetscReal block_width = (2.0 * PETSC_PI) / (PetscReal)blocks_per_dim;
107 PetscReal block_height = (2.0 * PETSC_PI) / (PetscReal)blocks_per_dim;
108
109 // Assign this block its specific sub-domain
110 user->Min_X = col * block_width;
111 user->Max_X = (col + 1) * block_width;
112 user->Min_Y = row * block_height;
113 user->Max_Y = (row + 1) * block_height;
114 user->Min_Z = 0.0;
115 user->Max_Z = 2.0 * PETSC_PI; // Z-domain is not decomposed
116 }
117 }
118 /*
119 * --- EXTENSIBILITY HOOK ---
120 * To add another analytical case with special grid requirements:
121 *
122 * else if (strcmp(simCtx->AnalyticalSolutionType, "ChannelFlow") == 0) {
123 * // ... implement logic to set domain for ChannelFlow case ...
124 * }
125 */
126 else {
127 // Fallback for other analytical solutions: require manual grid specification.
128 LOG(GLOBAL, LOG_INFO, "Analytical solution '%s' detected. Using user-provided grid generation inputs.\n", simCtx->AnalyticalSolutionType);
129 ierr = ReadGridGenerationInputs(user); CHKERRQ(ierr);
130 PetscFunctionReturn(0); // Return early as ReadGridGenerationInputs already handled everything.
131 }
132
133 // --- For TGV3D, we now read the grid resolution, as this is independent of the domain ---
134 PetscInt *IMs, *JMs, *KMs;
135 PetscBool found;
136 PetscInt count;
137
138 ierr = PetscMalloc3(nblk, &IMs, nblk, &JMs, nblk, &KMs); CHKERRQ(ierr);
139
140 // Set defaults first
141 for(PetscInt i=0; i<nblk; ++i) { IMs[i] = 10; JMs[i] = 10; KMs[i] = 10; }
142
143 count = nblk; ierr = PetscOptionsGetIntArray(NULL, NULL, "-im", IMs, &count, &found); CHKERRQ(ierr);
144 count = nblk; ierr = PetscOptionsGetIntArray(NULL, NULL, "-jm", JMs, &count, &found); CHKERRQ(ierr);
145 count = nblk; ierr = PetscOptionsGetIntArray(NULL, NULL, "-km", KMs, &count, &found); CHKERRQ(ierr);
146
147 user->IM = IMs[block_index];
148 user->JM = JMs[block_index];
149 user->KM = KMs[block_index];
150
151 // We can also read stretching ratios, as they are independent of the domain size
152 // For simplicity, we assume uniform grid unless specified.
153 user->rx = 1.0; user->ry = 1.0; user->rz = 1.0;
154
155 LOG_ALLOW(LOCAL, LOG_DEBUG, "Rank %d: Block %d grid resolution set: IM=%d, JM=%d, KM=%d\n",
156 simCtx->rank, block_index, user->IM, user->JM, user->KM);
157 LOG_ALLOW(LOCAL, LOG_DEBUG, "Rank %d: Block %d final bounds: X=[%.4f, %.4f], Y=[%.4f, %.4f], Z=[%.4f, %.4f]\n",
158 simCtx->rank, block_index, user->Min_X, user->Max_X, user->Min_Y, user->Max_Y, user->Min_Z, user->Max_Z);
159
160 ierr = PetscFree3(IMs, JMs, KMs); CHKERRQ(ierr);
161
163 PetscFunctionReturn(0);
164}
PetscErrorCode ReadGridGenerationInputs(UserCtx *user)
Parses command-line options for a programmatically generated grid for a SINGLE block.
Definition io.c:77
#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 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
#define PROFILE_FUNCTION_END
Marks the end of a profiled code block.
Definition logging.h:746
#define LOG(scope, level, fmt,...)
Logging macro for PETSc-based applications with scope control.
Definition logging.h:84
@ 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:588
PetscInt block_number
Definition variables.h:649
SimCtx * simCtx
Back-pointer to the master simulation context.
Definition variables.h:731
PetscReal Min_X
Definition variables.h:738
PetscInt KM
Definition variables.h:737
PetscInt _this
Definition variables.h:741
PetscReal ry
Definition variables.h:742
PetscReal Max_Y
Definition variables.h:738
PetscReal rz
Definition variables.h:742
PetscInt JM
Definition variables.h:737
PetscReal Min_Z
Definition variables.h:738
char AnalyticalSolutionType[PETSC_MAX_PATH_LEN]
Definition variables.h:617
PetscReal Max_X
Definition variables.h:738
PetscReal Min_Y
Definition variables.h:738
PetscInt IM
Definition variables.h:737
PetscReal rx
Definition variables.h:742
PetscReal Max_Z
Definition variables.h:738
The master context for the entire simulation.
Definition variables.h:585
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 from the simulation context and calls the corresponding private implementation function (e.g., for Taylor-Green Vortex, lid-driven cavity, etc.). This design keeps the main simulation code clean and makes it easy to add new analytical test cases.

Parameters
simCtxThe main simulation context, containing configuration and state.
Returns
PetscErrorCode 0 on success.

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 184 of file AnalyticalSolutions.c.

185{
186 PetscErrorCode ierr;
187 PetscFunctionBeginUser;
189
190 // -- Before any operation, here is defensive test to ensure that the Corner->Center Interpolation method works
191 //ierr = TestCornerToCenterInterpolation(&(simCtx->usermg.mgctx[simCtx->usermg.mglevels - 1]->user[0]));
192
193 // --- Dispatch based on the string provided by the user ---
194 if (strcmp(simCtx->AnalyticalSolutionType, "TGV3D") == 0) {
195 LOG_ALLOW(GLOBAL, LOG_DEBUG, "Applying Analytical Solution: 3D Taylor-Green Vortex (TGV3D).\n");
196 ierr = SetAnalyticalSolution_TGV3D(simCtx); CHKERRQ(ierr);
197 }
198 /*
199 * --- EXTENSIBILITY HOOK ---
200 * To add a new analytical solution (e.g., "ChannelFlow"):
201 * 1. Add an `else if` block here:
202 *
203 * else if (strcmp(simCtx->AnalyticalSolutionType, "ChannelFlow") == 0) {
204 * LOG_ALLOW(GLOBAL, LOG_DEBUG, "Applying Analytical Solution: Channel Flow.\n");
205 * ierr = SetAnalyticalSolution_ChannelFlow(simCtx); CHKERRQ(ierr);
206 * }
207 *
208 * 2. Implement the static function `SetAnalyticalSolution_ChannelFlow(SimCtx *simCtx)`
209 * below, following the TGV3D pattern.
210 */
211 else {
212 // If the type is unknown, raise a fatal error to prevent silent failures.
213 SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown AnalyticalSolutionType specified: '%s'", simCtx->AnalyticalSolutionType);
214 }
215
217 PetscFunctionReturn(0);
218}
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()

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 447 of file AnalyticalSolutions.c.

448{
449 PetscErrorCode ierr;
450 PetscInt nLocal;
451 PetscReal *vels;
452
453 PetscFunctionBeginUser;
454
455 LOG_ALLOW(GLOBAL, LOG_DEBUG, "Type: %s\n",
456 simCtx->AnalyticalSolutionType ? simCtx->AnalyticalSolutionType : "default");
457
458 // Check for specific analytical solution types
459 if (simCtx->AnalyticalSolutionType && strcmp(simCtx->AnalyticalSolutionType, "TGV3D") == 0) {
460 LOG_ALLOW(GLOBAL, LOG_DEBUG, "Using TGV3D solution.\n");
461 ierr = SetAnalyticalSolutionForParticles_TGV3D(tempVec, simCtx); CHKERRQ(ierr);
462 return 0;
463 }
464
465 ierr = VecRestoreArray(tempVec, &vels); CHKERRQ(ierr);
466
467 PetscFunctionReturn(0);
468}
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: