PICurv 0.1.0
A Parallel Particle-In-Cell Solver for Curvilinear LES
Loading...
Searching...
No Matches
AnalyticalSolutions.c
Go to the documentation of this file.
1/**
2 * @file AnalyticalSolutions.c
3 * @brief Implements the analytical solution engine for initializing or driving the simulation.
4 *
5 * @details This file provides a modular and extensible framework for applying analytical solutions
6 * to the Eulerian fields. The primary entry point is `AnalyticalSolutionEngine`, which acts
7 * as a dispatcher based on user configuration.
8 *
9 * --- DESIGN PHILOSOPHY ---
10 * 1. **Non-Dimensional Core:** All calculations within this engine are performed in
11 * **non-dimensional units** (e.g., reference velocity U_ref=1.0, reference length L_ref=1.0).
12 * This is critical for consistency with the core numerical solver, which also operates on
13 * non-dimensional equations. The `simCtx->scaling` parameters are intentionally NOT used here;
14 * they are reserved for dimensionalization during I/O and post-processing only.
15 *
16 * 2. **Separation of Concerns:** The role of this engine is to set the **physical state** of the
17 * fluid at a given time `t`. This involves:
18 * - For most types (TGV3D, ZERO_FLOW): setting `Ucat` and `P` directly.
19 * - For curvilinear-aware types (UNIFORM_FLOW): setting `Ucont` via metric
20 * dot products and deriving `Ucat` through `Contra2Cart`.
21 * - Declaring the physical values on the boundaries by populating the boundary condition
22 * vector (`user->Bcs.Ubcs`).
23 * It does NOT implement the numerical scheme for ghost cells. Instead, after setting the
24 * physical state, it relies on the solver's standard utility functions (`UpdateDummyCells`,
25 * `UpdateCornerNodes`) to correctly populate all ghost cell layers.
26 *
27 * 3. **Extensibility:** The dispatcher design makes it straightforward to add new analytical
28 * solutions. A developer only needs to add a new `else if` condition and a corresponding
29 * `static` implementation function, without modifying any other part of the solver.
30 */
31
32#include "AnalyticalSolutions.h" // The header that declares this file's functions
33#include <petscmath.h> // Provides PETSc-compatible math functions like sin, cos, exp
34
35// Forward-declare the private function for the TGV3D case.
36// This function is not visible outside this file, enforcing modularity.
37static PetscErrorCode SetAnalyticalSolution_TGV3D(SimCtx *simCtx);
38static PetscErrorCode SetAnalyticalSolution_ZeroFlow(SimCtx *simCtx);
39static PetscErrorCode SetAnalyticalSolution_UniformFlow(SimCtx *simCtx);
40static PetscErrorCode SetAnalyticalSolutionForParticles_TGV3D(Vec tempVec, SimCtx *simCtx);
41static PetscErrorCode SetAnalyticalSolutionForParticles_UniformFlow(Vec tempVec, SimCtx *simCtx);
42static PetscErrorCode EvaluateConfiguredScalarProfile(const VerificationScalarConfig *cfg,
43 PetscReal x,
44 PetscReal y,
45 PetscReal z,
46 PetscReal *value);
47
48/**
49 * @brief Implementation of \ref AnalyticalTypeRequiresCustomGeometry().
50 * @details Full API contract (arguments, ownership, side effects) is documented with
51 * the header declaration in `include/AnalyticalSolutions.h`.
52 * @see AnalyticalTypeRequiresCustomGeometry()
53 */
54
55PetscBool AnalyticalTypeRequiresCustomGeometry(const char *analytical_type)
56{
57 if (!analytical_type) return PETSC_FALSE;
58 return (strcmp(analytical_type, "TGV3D") == 0) ? PETSC_TRUE : PETSC_FALSE;
59}
60
61/**
62 * @brief Implementation of \ref AnalyticalTypeSupportsInterpolationError().
63 * @details Full API contract (arguments, ownership, side effects) is documented with
64 * the header declaration in `include/AnalyticalSolutions.h`.
65 * @see AnalyticalTypeSupportsInterpolationError()
66 */
67PetscBool AnalyticalTypeSupportsInterpolationError(const char *analytical_type)
68{
69 if (!analytical_type) return PETSC_FALSE;
70 if (strcmp(analytical_type, "ZERO_FLOW") == 0 || strcmp(analytical_type, "UNIFORM_FLOW") == 0) return PETSC_FALSE;
71 return PETSC_TRUE;
72}
73
74#undef __FUNCT__
75#define __FUNCT__ "SetAnalyticalGridInfo"
76/**
77 * @brief Internal helper implementation: `SetAnalyticalGridInfo()`.
78 * @details Local to this translation unit.
79 */
80PetscErrorCode SetAnalyticalGridInfo(UserCtx *user)
81{
82 SimCtx *simCtx = user->simCtx;
83 PetscInt nblk = simCtx->block_number;
84 PetscInt block_index = user->_this;
85
86 PetscFunctionBeginUser;
88
90 SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONGSTATE,
91 "SetAnalyticalGridInfo called for analytical type '%s' that does not require custom geometry.",
93 }
94 if (user->IM <= 0 || user->JM <= 0 || user->KM <= 0) {
95 SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONGSTATE,
96 "Analytical grid resolution is not initialized. Ensure IM/JM/KM are preloaded before SetAnalyticalGridInfo.");
97 }
98
99 if (strcmp(simCtx->AnalyticalSolutionType, "TGV3D") == 0) {
100 LOG_ALLOW_SYNC(GLOBAL, LOG_DEBUG, "Rank %d: Configuring grid for TGV3D analytical solution, block %d.\n", simCtx->rank, block_index);
101
102 if (nblk == 1) {
103 // --- Single Block Case ---
104 if (block_index == 0) {
105 LOG_ALLOW(GLOBAL, LOG_INFO, "Single block detected. Setting domain to [0, 2*PI].\n");
106 }
107 user->Min_X = 0.0; user->Max_X = 2.0 * PETSC_PI;
108 user->Min_Y = 0.0; user->Max_Y = 2.0 * PETSC_PI;
109 user->Min_Z = 0.0; user->Max_Z = 0.2 * PETSC_PI; //2.0 * PETSC_PI;
110
111 } else { // --- Multi-Block Case ---
112 PetscReal s = sqrt((PetscReal)nblk);
113
114 // Validate that nblk is a perfect square.
115 if (fabs(s - floor(s)) > 1e-9) {
116 SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_INCOMP,
117 "\n\n*** CONFIGURATION ERROR FOR TGV3D ***\n"
118 "For multi-block TGV3D cases, the number of blocks must be a perfect square (e.g., 4, 9, 16).\n"
119 "You have specified %d blocks. Please adjust `-block_number`.\n", nblk);
120 }
121 PetscInt blocks_per_dim = (PetscInt)s;
122
123 if (block_index == 0) {
124 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);
125 }
126
127 // Determine the (row, col) position of this block in the 2D decomposition
128 PetscInt row = block_index / blocks_per_dim;
129 PetscInt col = block_index % blocks_per_dim;
130
131 // Calculate the width/height of each sub-domain
132 PetscReal block_width = (2.0 * PETSC_PI) / (PetscReal)blocks_per_dim;
133 PetscReal block_height = (2.0 * PETSC_PI) / (PetscReal)blocks_per_dim;
134
135 // Assign this block its specific sub-domain
136 user->Min_X = col * block_width;
137 user->Max_X = (col + 1) * block_width;
138 user->Min_Y = row * block_height;
139 user->Max_Y = (row + 1) * block_height;
140 user->Min_Z = 0.0;
141 user->Max_Z = 2.0 * PETSC_PI; // Z-domain is not decomposed
142 }
143 }
144 /*
145 * --- EXTENSIBILITY HOOK ---
146 * To add another analytical case with special grid requirements:
147 *
148 * else if (strcmp(simCtx->AnalyticalSolutionType, "ChannelFlow") == 0) {
149 * // ... implement logic to set domain for ChannelFlow case ...
150 * }
151 */
152 else {
153 SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_UNKNOWN_TYPE,
154 "Analytical type '%s' has no custom geometry implementation.",
155 simCtx->AnalyticalSolutionType);
156 }
157
158 // We can also read stretching ratios, as they are independent of the domain size
159 // For simplicity, we assume uniform grid unless specified.
160 user->rx = 1.0; user->ry = 1.0; user->rz = 1.0;
161
162 LOG_ALLOW(LOCAL, LOG_DEBUG, "Rank %d: Block %d grid resolution set: IM=%d, JM=%d, KM=%d\n",
163 simCtx->rank, block_index, user->IM, user->JM, user->KM);
164 LOG_ALLOW(LOCAL, LOG_DEBUG, "Rank %d: Block %d final bounds: X=[%.4f, %.4f], Y=[%.4f, %.4f], Z=[%.4f, %.4f]\n",
165 simCtx->rank, block_index, user->Min_X, user->Max_X, user->Min_Y, user->Max_Y, user->Min_Z, user->Max_Z);
166
168 PetscFunctionReturn(0);
169}
170
171
172/*================================================================================*
173 * PUBLIC ANALYTICAL SOLUTION ENGINE *
174 *================================================================================*/
175
176#undef __FUNCT__
177#define __FUNCT__ "AnalyticalSolutionEngine"
178/**
179 * @brief Implementation of \ref AnalyticalSolutionEngine().
180 * @details Full API contract (arguments, ownership, side effects) is documented with
181 * the header declaration in `include/AnalyticalSolutions.h`.
182 * @see AnalyticalSolutionEngine()
183 */
184PetscErrorCode AnalyticalSolutionEngine(SimCtx *simCtx)
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 else if (strcmp(simCtx->AnalyticalSolutionType, "ZERO_FLOW") == 0) {
199 LOG_ALLOW(GLOBAL, LOG_DEBUG, "Applying Analytical Solution: Zero Background Flow (ZERO_FLOW).\n");
200 ierr = SetAnalyticalSolution_ZeroFlow(simCtx); CHKERRQ(ierr);
201 }
202 else if (strcmp(simCtx->AnalyticalSolutionType, "UNIFORM_FLOW") == 0) {
203 LOG_ALLOW(GLOBAL, LOG_DEBUG, "Applying Analytical Solution: Uniform Background Flow (UNIFORM_FLOW).\n");
204 ierr = SetAnalyticalSolution_UniformFlow(simCtx); CHKERRQ(ierr);
205 }
206 /*
207 * --- EXTENSIBILITY HOOK ---
208 * To add a new analytical solution (e.g., "ChannelFlow"):
209 * 1. Add an `else if` block here:
210 *
211 * else if (strcmp(simCtx->AnalyticalSolutionType, "ChannelFlow") == 0) {
212 * LOG_ALLOW(GLOBAL, LOG_DEBUG, "Applying Analytical Solution: Channel Flow.\n");
213 * ierr = SetAnalyticalSolution_ChannelFlow(simCtx); CHKERRQ(ierr);
214 * }
215 *
216 * 2. Implement the static function `SetAnalyticalSolution_ChannelFlow(SimCtx *simCtx)`
217 * below, following the TGV3D pattern.
218 */
219 else {
220 // If the type is unknown, raise a fatal error to prevent silent failures.
221 SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown AnalyticalSolutionType specified: '%s'", simCtx->AnalyticalSolutionType);
222 }
223
225 PetscFunctionReturn(0);
226}
227
228
229/*================================================================================*
230 * PRIVATE IMPLEMENTATIONS *
231 *================================================================================*/
232
233#undef __FUNCT__
234#define __FUNCT__ "SetAnalyticalSolution_TGV3D"
235/**
236 * @brief Internal helper implementation: `SetAnalyticalSolution_TGV3D()`.
237 * @details Local to this translation unit.
238 */
239static PetscErrorCode SetAnalyticalSolution_TGV3D(SimCtx *simCtx)
240{
241 PetscErrorCode ierr;
242 UserCtx *user_finest = simCtx->usermg.mgctx[simCtx->usermg.mglevels - 1].user;
243
244 // --- NON-DIMENSIONAL TGV Parameters ---
245 const PetscReal V0 = 1.0; // Non-dimensional reference velocity.
246 const PetscReal rho = 1.0; // Non-dimensional reference density.
247 const PetscReal p0 = 0.0; // Non-dimensional reference pressure.
248
249 // Kinematic viscosity is derived from the non-dimensional Reynolds number.
250 const PetscReal nu = (simCtx->ren > 0) ? (1.0 / simCtx->ren) : 0.0;
251
252 const PetscReal k = 1.0; // Wavenumber, assumes a non-dimensional [0, 2*pi] domain.
253 const PetscReal t = simCtx->ti;
254
255 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);
256
257 const PetscReal vel_decay = exp(-2.0 * nu * k * k * t);
258 const PetscReal prs_decay = exp(-4.0 * nu * k * k * t);
259
260 PetscFunctionBeginUser;
261
262 for (PetscInt bi = 0; bi < simCtx->block_number; bi++) {
263 UserCtx* user = &user_finest[bi];
264 DMDALocalInfo info = user->info;
265 PetscInt xs = info.xs, xe = info.xs + info.xm;
266 PetscInt ys = info.ys, ye = info.ys + info.ym;
267 PetscInt zs = info.zs, ze = info.zs + info.zm;
268 PetscInt mx = info.mx, my = info.my, mz = info.mz;
269
270 Cmpnts ***ucat, ***ubcs;
271 const Cmpnts ***cent, ***cent_x, ***cent_y, ***cent_z;
272 PetscReal ***p;
273
274 // Define loop bounds for physical interior cells owned by this rank.
275 PetscInt lxs = (xs == 0) ? xs + 1 : xs, lxe = (xe == mx) ? xe - 1 : xe;
276 PetscInt lys = (ys == 0) ? ys + 1 : ys, lye = (ye == my) ? ye - 1 : ye;
277 PetscInt lzs = (zs == 0) ? zs + 1 : zs, lze = (ze == mz) ? ze - 1 : ze;
278
279 // --- Get Arrays ---
280 ierr = DMDAVecGetArray(user->fda, user->Ucat, &ucat); CHKERRQ(ierr);
281 ierr = DMDAVecGetArray(user->da, user->P, &p); CHKERRQ(ierr);
282 ierr = DMDAVecGetArray(user->fda, user->Bcs.Ubcs, &ubcs); CHKERRQ(ierr);
283 ierr = DMDAVecGetArrayRead(user->fda, user->Cent, &cent); CHKERRQ(ierr);
284 ierr = DMDAVecGetArrayRead(user->fda, user->Centx, &cent_x); CHKERRQ(ierr);
285 ierr = DMDAVecGetArrayRead(user->fda, user->Centy, &cent_y); CHKERRQ(ierr);
286 ierr = DMDAVecGetArrayRead(user->fda, user->Centz, &cent_z); CHKERRQ(ierr);
287
288 // --- Set INTERIOR cell-centered velocity (Ucat) ---
289 for (PetscInt k_cell = lzs; k_cell < lze; k_cell++) {
290 for (PetscInt j_cell = lys; j_cell < lye; j_cell++) {
291 for (PetscInt i_cell = lxs; i_cell < lxe; i_cell++) {
292 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;
293 ucat[k_cell][j_cell][i_cell].x = V0 * sin(k*cx) * cos(k*cy) * cos(k*cz) * vel_decay;
294 ucat[k_cell][j_cell][i_cell].y = -V0 * cos(k*cx) * sin(k*cy) * cos(k*cz) * vel_decay;
295 ucat[k_cell][j_cell][i_cell].z = 0.0;
296 }
297 }
298 }
299
300
301 // --- Set INTERIOR cell-centered pressure (P) ---
302 for (PetscInt k_cell = lzs; k_cell < lze; k_cell++) {
303 for (PetscInt j_cell = lys; j_cell < lye; j_cell++) {
304 for (PetscInt i_cell = lxs; i_cell < lxe; i_cell++) {
305 const PetscReal cx = cent[k_cell][j_cell][i_cell].x, cy = cent[k_cell][j_cell][i_cell].y;
306 p[k_cell][j_cell][i_cell] = p0 + (rho * V0 * V0 / 4.0) * (cos(2*k*cx) + cos(2*k*cy)) * prs_decay;
307 }
308 }
309 }
310
311
312 // --- Set BOUNDARY condition vector for velocity (Ubcs) ---
313 if (xs == 0) for (PetscInt k=zs; k<ze; k++) for (PetscInt j=ys; j<ye; j++) {
314 const PetscReal fcx=cent_x[k][j][xs].x, fcy=cent_x[k][j][xs].y, fcz=cent_x[k][j][xs].z;
315 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;
316 }
317 if (xe == mx) for (PetscInt k=zs; k<ze; k++) for (PetscInt j=ys; j<ye; j++) {
318 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;
319 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;
320 }
321 if (ys == 0) for (PetscInt k=zs; k<ze; k++) for (PetscInt i=xs; i<xe; i++) {
322 const PetscReal fcx=cent_y[k][ys][i].x, fcy=cent_y[k][ys][i].y, fcz=cent_y[k][ys][i].z;
323 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;
324 }
325 if (ye == my) for (PetscInt k=zs; k<ze; k++) for (PetscInt i=xs; i<xe; i++) {
326 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;
327 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;
328 }
329 if (zs == 0) for (PetscInt j=ys; j<ye; j++) for (PetscInt i=xs; i<xe; i++) {
330 const PetscReal fcx=cent_z[zs][j][i].x, fcy=cent_z[zs][j][i].y, fcz=cent_z[zs][j][i].z;
331 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;
332 }
333 if (ze == mz) for (PetscInt j=ys; j<ye; j++) for (PetscInt i=xs; i<xe; i++) {
334 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;
335 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;
336 }
337
338 // --- Set PRESSURE GHOST CELLS (Neumann BC: P_ghost = P_interior) ---
339 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];
340 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];
341
342 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];
343 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];
344
345 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];
346 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];
347
348 // --- Restore all arrays ---
349 ierr = DMDAVecRestoreArray(user->fda, user->Ucat, &ucat); CHKERRQ(ierr);
350 ierr = DMDAVecRestoreArray(user->da, user->P, &p); CHKERRQ(ierr);
351 ierr = DMDAVecRestoreArray(user->fda, user->Bcs.Ubcs, &ubcs); CHKERRQ(ierr);
352 ierr = DMDAVecRestoreArrayRead(user->fda, user->Cent, &cent); CHKERRQ(ierr);
353 ierr = DMDAVecRestoreArrayRead(user->fda, user->Centx, &cent_x); CHKERRQ(ierr);
354 ierr = DMDAVecRestoreArrayRead(user->fda, user->Centy, &cent_y); CHKERRQ(ierr);
355 ierr = DMDAVecRestoreArrayRead(user->fda, user->Centz, &cent_z); CHKERRQ(ierr);
356
357 // Pre-Dummy cell update synchronization.
358 ierr = UpdateLocalGhosts(user,"Ucat");
359 ierr = UpdateLocalGhosts(user,"P");
360
361 // --- Finalize all ghost cell values ---
362 ierr = UpdateDummyCells(user); CHKERRQ(ierr);
363 ierr = UpdateCornerNodes(user); CHKERRQ(ierr);
364
365 // Final Synchronization.
366 ierr = UpdateLocalGhosts(user,"Ucat");
367 ierr = UpdateLocalGhosts(user,"P");
368
369 }
370
371 PetscFunctionReturn(0);
372}
373
374#undef __FUNCT__
375#define __FUNCT__ "SetAnalyticalSolution_ZeroFlow"
376/**
377 * @brief Internal helper implementation: `SetAnalyticalSolution_ZeroFlow()`.
378 * @details Local to this translation unit.
379 */
380static PetscErrorCode SetAnalyticalSolution_ZeroFlow(SimCtx *simCtx)
381{
382 PetscErrorCode ierr;
383 UserCtx *user_finest = simCtx->usermg.mgctx[simCtx->usermg.mglevels - 1].user;
384
385 PetscFunctionBeginUser;
386
387 for (PetscInt bi = 0; bi < simCtx->block_number; bi++) {
388 UserCtx *user = &user_finest[bi];
389
390 ierr = VecZeroEntries(user->Ucat); CHKERRQ(ierr);
391 ierr = VecZeroEntries(user->P); CHKERRQ(ierr);
392 ierr = VecZeroEntries(user->Bcs.Ubcs); CHKERRQ(ierr);
393
394 // Ghost-cell finalization — identical sequence to TGV3D
395 ierr = UpdateLocalGhosts(user, "Ucat"); CHKERRQ(ierr);
396 ierr = UpdateLocalGhosts(user, "P"); CHKERRQ(ierr);
397 ierr = UpdateDummyCells(user); CHKERRQ(ierr);
398 ierr = UpdateCornerNodes(user); CHKERRQ(ierr);
399 ierr = UpdateLocalGhosts(user, "Ucat"); CHKERRQ(ierr);
400 ierr = UpdateLocalGhosts(user, "P"); CHKERRQ(ierr);
401 }
402
403 PetscFunctionReturn(0);
404}
405
406#undef __FUNCT__
407#define __FUNCT__ "SetAnalyticalSolution_UniformFlow"
408/**
409 * @brief Internal helper implementation: `SetAnalyticalSolution_UniformFlow()`.
410 * @details Local to this translation unit.
411 */
412static PetscErrorCode SetAnalyticalSolution_UniformFlow(SimCtx *simCtx)
413{
414 PetscErrorCode ierr;
415 UserCtx *user_finest = simCtx->usermg.mgctx[simCtx->usermg.mglevels - 1].user;
416 const Cmpnts uniform_velocity = simCtx->AnalyticalUniformVelocity;
417 const PetscReal u = uniform_velocity.x;
418 const PetscReal v = uniform_velocity.y;
419 const PetscReal w = uniform_velocity.z;
420
421 PetscFunctionBeginUser;
422
423 for (PetscInt bi = 0; bi < simCtx->block_number; bi++) {
424 UserCtx *user = &user_finest[bi];
425 DMDALocalInfo info = user->info;
426 PetscInt xs = info.xs, xe = info.xs + info.xm;
427 PetscInt ys = info.ys, ye = info.ys + info.ym;
428 PetscInt zs = info.zs, ze = info.zs + info.zm;
429 PetscInt mx = info.mx, my = info.my, mz = info.mz;
430
431 // --- Step 1: Set Ucont (contravariant flux) using metric dot products ---
432 // U^i = g_i · u_phys, where g_i is the face area vector (Csi, Eta, Zet).
433 // This is the curvilinear form: the volume flux through each face.
434 // We loop over ALL owned nodes (including boundary faces) because Contra2Cart
435 // averages adjacent face values — boundary faces must carry valid fluxes.
436 Cmpnts ***ucont, ***ubcs;
437 const Cmpnts ***csi, ***eta, ***zet;
438
439 ierr = DMDAVecGetArray(user->fda, user->Ucont, &ucont); CHKERRQ(ierr);
440 ierr = DMDAVecGetArray(user->fda, user->Bcs.Ubcs, &ubcs); CHKERRQ(ierr);
441 ierr = DMDAVecGetArrayRead(user->fda, user->lCsi, &csi); CHKERRQ(ierr);
442 ierr = DMDAVecGetArrayRead(user->fda, user->lEta, &eta); CHKERRQ(ierr);
443 ierr = DMDAVecGetArrayRead(user->fda, user->lZet, &zet); CHKERRQ(ierr);
444
445 for (PetscInt k = zs; k < ze; k++) {
446 for (PetscInt j = ys; j < ye; j++) {
447 for (PetscInt i = xs; i < xe; i++) {
448 ucont[k][j][i].x = csi[k][j][i].x * u + csi[k][j][i].y * v + csi[k][j][i].z * w;
449 ucont[k][j][i].y = eta[k][j][i].x * u + eta[k][j][i].y * v + eta[k][j][i].z * w;
450 ucont[k][j][i].z = zet[k][j][i].x * u + zet[k][j][i].y * v + zet[k][j][i].z * w;
451 }
452 }
453 }
454
455 // --- Step 2: Set Ubcs at boundaries (physical Cartesian velocity) ---
456 if (xs == 0) for (PetscInt k = zs; k < ze; k++) for (PetscInt j = ys; j < ye; j++) ubcs[k][j][xs] = uniform_velocity;
457 if (xe == mx) for (PetscInt k = zs; k < ze; k++) for (PetscInt j = ys; j < ye; j++) ubcs[k][j][xe - 1] = uniform_velocity;
458 if (ys == 0) for (PetscInt k = zs; k < ze; k++) for (PetscInt i = xs; i < xe; i++) ubcs[k][ys][i] = uniform_velocity;
459 if (ye == my) for (PetscInt k = zs; k < ze; k++) for (PetscInt i = xs; i < xe; i++) ubcs[k][ye - 1][i] = uniform_velocity;
460 if (zs == 0) for (PetscInt j = ys; j < ye; j++) for (PetscInt i = xs; i < xe; i++) ubcs[zs][j][i] = uniform_velocity;
461 if (ze == mz) for (PetscInt j = ys; j < ye; j++) for (PetscInt i = xs; i < xe; i++) ubcs[ze - 1][j][i] = uniform_velocity;
462
463 ierr = DMDAVecRestoreArrayRead(user->fda, user->lZet, &zet); CHKERRQ(ierr);
464 ierr = DMDAVecRestoreArrayRead(user->fda, user->lEta, &eta); CHKERRQ(ierr);
465 ierr = DMDAVecRestoreArrayRead(user->fda, user->lCsi, &csi); CHKERRQ(ierr);
466 ierr = DMDAVecRestoreArray(user->fda, user->Bcs.Ubcs, &ubcs); CHKERRQ(ierr);
467 ierr = DMDAVecRestoreArray(user->fda, user->Ucont, &ucont); CHKERRQ(ierr);
468
469 // --- Step 3: Zero pressure ---
470 ierr = VecZeroEntries(user->P); CHKERRQ(ierr);
471
472 // --- Step 4: Finalize state — derive Ucat from Ucont via metric inversion ---
473 ierr = UpdateLocalGhosts(user, "Ucont"); CHKERRQ(ierr);
474 ierr = Contra2Cart(user); CHKERRQ(ierr);
475 ierr = UpdateLocalGhosts(user, "Ucat"); CHKERRQ(ierr);
476 ierr = UpdateLocalGhosts(user, "P"); CHKERRQ(ierr);
477 ierr = UpdateDummyCells(user); CHKERRQ(ierr);
478 ierr = UpdateCornerNodes(user); CHKERRQ(ierr);
479 ierr = UpdateLocalGhosts(user, "Ucat"); CHKERRQ(ierr);
480 ierr = UpdateLocalGhosts(user, "P"); CHKERRQ(ierr);
481 }
482
483 PetscFunctionReturn(0);
484}
485
486#undef __FUNCT__
487#define __FUNCT__ "SetAnalyticalSolutionForParticles_TGV3D"
488/**
489 * @brief Internal helper implementation: `SetAnalyticalSolutionForParticles_TGV3D()`.
490 * @details Local to this translation unit.
491 */
492static PetscErrorCode SetAnalyticalSolutionForParticles_TGV3D(Vec tempVec, SimCtx *simCtx)
493{
494 PetscErrorCode ierr;
495 PetscInt nLocal;
496 PetscReal *data;
497
498 PetscFunctionBeginUser;
499
500 // TGV3D parameters (matching your Eulerian implementation)
501 const PetscReal V0 = 1.0;
502 const PetscReal k = 1.0;
503 const PetscReal nu = (simCtx->ren > 0) ? (1.0 / simCtx->ren) : 0.0;
504 const PetscReal t = simCtx->ti;
505 const PetscReal vel_decay = exp(-2.0 * nu * k * k * t);
506
507 LOG_ALLOW(GLOBAL, LOG_DEBUG, "TGV3D Particles: t=%.4f, V0=%.4f, k=%.4f, nu=%.6f\n", t, V0, k, nu);
508
509 ierr = VecGetLocalSize(tempVec, &nLocal); CHKERRQ(ierr);
510 ierr = VecGetArray(tempVec, &data); CHKERRQ(ierr);
511
512 // Process particles: data is interleaved [x0,y0,z0, x1,y1,z1, ...]
513 for (PetscInt i = 0; i < nLocal; i += 3) {
514 const PetscReal x = data[i];
515 const PetscReal y = data[i+1];
516 const PetscReal z = data[i+2];
517
518 // TGV3D velocity field
519 data[i] = V0 * sin(k*x) * cos(k*y) * cos(k*z) * vel_decay; // u
520 data[i+1] = -V0 * cos(k*x) * sin(k*y) * cos(k*z) * vel_decay; // v
521 data[i+2] = 0.0; // w
522 }
523
524 ierr = VecRestoreArray(tempVec, &data); CHKERRQ(ierr);
525
526 PetscFunctionReturn(0);
527}
528
529#undef __FUNCT__
530#define __FUNCT__ "SetAnalyticalSolutionForParticles_UniformFlow"
531/**
532 * @brief Internal helper implementation: `SetAnalyticalSolutionForParticles_UniformFlow()`.
533 * @details Local to this translation unit.
534 */
535static PetscErrorCode SetAnalyticalSolutionForParticles_UniformFlow(Vec tempVec, SimCtx *simCtx)
536{
537 PetscErrorCode ierr;
538 PetscInt nLocal;
539 PetscReal *data;
540 const Cmpnts uniform_velocity = simCtx->AnalyticalUniformVelocity;
541
542 PetscFunctionBeginUser;
543
544 ierr = VecGetLocalSize(tempVec, &nLocal); CHKERRQ(ierr);
545 ierr = VecGetArray(tempVec, &data); CHKERRQ(ierr);
546 for (PetscInt i = 0; i < nLocal; i += 3) {
547 data[i] = uniform_velocity.x;
548 data[i + 1] = uniform_velocity.y;
549 data[i + 2] = uniform_velocity.z;
550 }
551 ierr = VecRestoreArray(tempVec, &data); CHKERRQ(ierr);
552
553 PetscFunctionReturn(0);
554}
555
556#undef __FUNCT__
557#define __FUNCT__ "SetAnalyticalSolutionForParticles"
558/**
559 * @brief Implementation of \ref SetAnalyticalSolutionForParticles().
560 * @details Full API contract (arguments, ownership, side effects) is documented with
561 * the header declaration in `include/AnalyticalSolutions.h`.
562 * @see SetAnalyticalSolutionForParticles()
563 */
564PetscErrorCode SetAnalyticalSolutionForParticles(Vec tempVec, SimCtx *simCtx)
565{
566 PetscErrorCode ierr;
567 const char *analytical_type = simCtx->AnalyticalSolutionType[0] ? simCtx->AnalyticalSolutionType : "default";
568
569 PetscFunctionBeginUser;
570
571 LOG_ALLOW(GLOBAL, LOG_DEBUG, "Type: %s\n", analytical_type);
572
573 // Check for specific analytical solution types
574 if (strcmp(simCtx->AnalyticalSolutionType, "TGV3D") == 0) {
575 LOG_ALLOW(GLOBAL, LOG_DEBUG, "Using TGV3D solution.\n");
576 ierr = SetAnalyticalSolutionForParticles_TGV3D(tempVec, simCtx); CHKERRQ(ierr);
577 return 0;
578 }
579 if (strcmp(simCtx->AnalyticalSolutionType, "UNIFORM_FLOW") == 0) {
580 LOG_ALLOW(GLOBAL, LOG_DEBUG, "Using UNIFORM_FLOW solution.\n");
581 ierr = SetAnalyticalSolutionForParticles_UniformFlow(tempVec, simCtx); CHKERRQ(ierr);
582 return 0;
583 }
584
585 PetscFunctionReturn(0);
586}
587
588/**
589 * @brief Internal helper that evaluates the configured scalar verification profile.
590 * @details Local to this translation unit.
591 */
593 PetscReal x,
594 PetscReal y,
595 PetscReal z,
596 PetscReal *value)
597{
598 PetscFunctionBeginUser;
599 if (!cfg) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "VerificationScalarConfig cannot be NULL.");
600 if (!value) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "Scalar output pointer cannot be NULL.");
601
602 if (strcmp(cfg->profile, "CONSTANT") == 0) {
603 *value = cfg->value;
604 } else if (strcmp(cfg->profile, "LINEAR_X") == 0) {
605 *value = cfg->phi0 + cfg->slope_x * x;
606 } else if (strcmp(cfg->profile, "SIN_PRODUCT") == 0) {
607 *value = cfg->amplitude *
608 PetscSinReal(cfg->kx * x) *
609 PetscSinReal(cfg->ky * y) *
610 PetscSinReal(cfg->kz * z);
611 } else {
612 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG,
613 "Unsupported verification scalar profile '%s'.", cfg->profile);
614 }
615
616 PetscFunctionReturn(0);
617}
618
619#undef __FUNCT__
620#define __FUNCT__ "EvaluateAnalyticalScalarProfile"
621/**
622 * @brief Implementation of \ref EvaluateAnalyticalScalarProfile().
623 * @details Full API contract (arguments, ownership, side effects) is documented with
624 * the header declaration in `include/AnalyticalSolutions.h`.
625 * @see EvaluateAnalyticalScalarProfile()
626 */
627PetscErrorCode EvaluateAnalyticalScalarProfile(const SimCtx *simCtx,
628 PetscReal x,
629 PetscReal y,
630 PetscReal z,
631 PetscReal t,
632 PetscReal *value)
633{
634 PetscFunctionBeginUser;
635 (void)t;
636 if (!simCtx) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "SimCtx cannot be NULL.");
637 if (!simCtx->verificationScalar.enabled) {
638 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE,
639 "EvaluateAnalyticalScalarProfile requires verification scalar mode to be enabled.");
640 }
641 PetscCall(EvaluateConfiguredScalarProfile(&simCtx->verificationScalar, x, y, z, value));
642 PetscFunctionReturn(0);
643}
644
645#undef __FUNCT__
646#define __FUNCT__ "SetAnalyticalScalarFieldOnParticles"
647/**
648 * @brief Implementation of \ref SetAnalyticalScalarFieldOnParticles().
649 * @details Full API contract (arguments, ownership, side effects) is documented with
650 * the header declaration in `include/AnalyticalSolutions.h`.
651 * @see SetAnalyticalScalarFieldOnParticles()
652 */
653PetscErrorCode SetAnalyticalScalarFieldOnParticles(UserCtx *user, const char *swarm_field_name)
654{
655 PetscErrorCode ierr;
656 PetscInt nlocal = 0;
657 PetscReal *positions = NULL;
658 PetscReal *scalar_values = NULL;
659
660 PetscFunctionBeginUser;
661 if (!user) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "UserCtx cannot be NULL.");
662 if (!user->swarm) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "UserCtx->swarm is NULL.");
663 if (!swarm_field_name) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "swarm_field_name cannot be NULL.");
664
665 ierr = DMSwarmGetLocalSize(user->swarm, &nlocal); CHKERRQ(ierr);
666 if (nlocal == 0) PetscFunctionReturn(0);
667
668 ierr = DMSwarmGetField(user->swarm, "position", NULL, NULL, (void **)&positions); CHKERRQ(ierr);
669 ierr = DMSwarmGetField(user->swarm, swarm_field_name, NULL, NULL, (void **)&scalar_values); CHKERRQ(ierr);
670
671 for (PetscInt p = 0; p < nlocal; ++p) {
672 PetscReal value = 0.0;
674 positions[3 * p + 0],
675 positions[3 * p + 1],
676 positions[3 * p + 2],
677 user->simCtx->ti,
678 &value); CHKERRQ(ierr);
679 scalar_values[p] = value;
680 }
681
682 ierr = DMSwarmRestoreField(user->swarm, swarm_field_name, NULL, NULL, (void **)&scalar_values); CHKERRQ(ierr);
683 ierr = DMSwarmRestoreField(user->swarm, "position", NULL, NULL, (void **)&positions); CHKERRQ(ierr);
684 PetscFunctionReturn(0);
685}
686
687#undef __FUNCT__
688#define __FUNCT__ "SetAnalyticalScalarFieldAtCellCenters"
689/**
690 * @brief Implementation of \ref SetAnalyticalScalarFieldAtCellCenters().
691 * @details Full API contract (arguments, ownership, side effects) is documented with
692 * the header declaration in `include/AnalyticalSolutions.h`.
693 * @see SetAnalyticalScalarFieldAtCellCenters()
694 */
695PetscErrorCode SetAnalyticalScalarFieldAtCellCenters(UserCtx *user, Vec targetVec)
696{
697 PetscErrorCode ierr;
698 PetscReal ***target = NULL;
699 const Cmpnts ***cent = NULL;
700 DMDALocalInfo info;
701 PetscInt xs, xe, ys, ye, zs, ze, mx, my, mz;
702 PetscInt lxs, lxe, lys, lye, lzs, lze;
703
704 PetscFunctionBeginUser;
705 if (!user) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "UserCtx cannot be NULL.");
706 if (!targetVec) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "targetVec cannot be NULL.");
707
708 info = user->info;
709 xs = info.xs; xe = info.xs + info.xm;
710 ys = info.ys; ye = info.ys + info.ym;
711 zs = info.zs; ze = info.zs + info.zm;
712 mx = info.mx; my = info.my; mz = info.mz;
713 lxs = (xs == 0) ? xs + 1 : xs; lxe = (xe == mx) ? xe - 1 : xe;
714 lys = (ys == 0) ? ys + 1 : ys; lye = (ye == my) ? ye - 1 : ye;
715 lzs = (zs == 0) ? zs + 1 : zs; lze = (ze == mz) ? ze - 1 : ze;
716
717 ierr = VecSet(targetVec, 0.0); CHKERRQ(ierr);
718 ierr = DMDAVecGetArray(user->da, targetVec, &target); CHKERRQ(ierr);
719 ierr = DMDAVecGetArrayRead(user->fda, user->Cent, &cent); CHKERRQ(ierr);
720
721 for (PetscInt k = lzs; k < lze; ++k) {
722 for (PetscInt j = lys; j < lye; ++j) {
723 for (PetscInt i = lxs; i < lxe; ++i) {
724 PetscReal value = 0.0;
726 cent[k][j][i].x,
727 cent[k][j][i].y,
728 cent[k][j][i].z,
729 user->simCtx->ti,
730 &value); CHKERRQ(ierr);
731 target[k][j][i] = value;
732 }
733 }
734 }
735
736 ierr = DMDAVecRestoreArrayRead(user->fda, user->Cent, &cent); CHKERRQ(ierr);
737 ierr = DMDAVecRestoreArray(user->da, targetVec, &target); CHKERRQ(ierr);
738 PetscFunctionReturn(0);
739}
static PetscErrorCode SetAnalyticalSolution_ZeroFlow(SimCtx *simCtx)
Internal helper implementation: SetAnalyticalSolution_ZeroFlow().
PetscErrorCode SetAnalyticalScalarFieldOnParticles(UserCtx *user, const char *swarm_field_name)
Implementation of SetAnalyticalScalarFieldOnParticles().
static PetscErrorCode SetAnalyticalSolution_UniformFlow(SimCtx *simCtx)
Internal helper implementation: SetAnalyticalSolution_UniformFlow().
PetscErrorCode EvaluateAnalyticalScalarProfile(const SimCtx *simCtx, PetscReal x, PetscReal y, PetscReal z, PetscReal t, PetscReal *value)
Implementation of EvaluateAnalyticalScalarProfile().
PetscErrorCode SetAnalyticalScalarFieldAtCellCenters(UserCtx *user, Vec targetVec)
Implementation of SetAnalyticalScalarFieldAtCellCenters().
static PetscErrorCode SetAnalyticalSolutionForParticles_UniformFlow(Vec tempVec, SimCtx *simCtx)
Internal helper implementation: SetAnalyticalSolutionForParticles_UniformFlow().
static PetscErrorCode EvaluateConfiguredScalarProfile(const VerificationScalarConfig *cfg, PetscReal x, PetscReal y, PetscReal z, PetscReal *value)
Internal helper that evaluates the configured scalar verification profile.
PetscErrorCode AnalyticalSolutionEngine(SimCtx *simCtx)
Implementation of AnalyticalSolutionEngine().
PetscBool AnalyticalTypeRequiresCustomGeometry(const char *analytical_type)
Implementation of AnalyticalTypeRequiresCustomGeometry().
PetscBool AnalyticalTypeSupportsInterpolationError(const char *analytical_type)
Implementation of AnalyticalTypeSupportsInterpolationError().
PetscErrorCode SetAnalyticalSolutionForParticles(Vec tempVec, SimCtx *simCtx)
Implementation of SetAnalyticalSolutionForParticles().
static PetscErrorCode SetAnalyticalSolution_TGV3D(SimCtx *simCtx)
Internal helper implementation: SetAnalyticalSolution_TGV3D().
PetscErrorCode SetAnalyticalGridInfo(UserCtx *user)
Internal helper implementation: SetAnalyticalGridInfo().
static PetscErrorCode SetAnalyticalSolutionForParticles_TGV3D(Vec tempVec, SimCtx *simCtx)
Internal helper implementation: SetAnalyticalSolutionForParticles_TGV3D().
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 LOG_ALLOW_SYNC(scope, level, fmt,...)
Synchronized logging macro that checks both the log level and whether the calling function is in the ...
Definition logging.h:252
#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:199
#define PROFILE_FUNCTION_END
Marks the end of a profiled code block.
Definition logging.h:779
@ LOG_TRACE
Very fine-grained tracing information for in-depth debugging.
Definition logging.h:32
@ LOG_INFO
Informational messages about program execution.
Definition logging.h:30
@ LOG_DEBUG
Detailed debugging information.
Definition logging.h:31
#define PROFILE_FUNCTION_BEGIN
Marks the beginning of a profiled code block (typically a function).
Definition logging.h:770
PetscErrorCode Contra2Cart(UserCtx *user)
Reconstructs Cartesian velocity (Ucat) at cell centers from contravariant velocity (Ucont) defined on...
Definition setup.c:2247
PetscErrorCode UpdateLocalGhosts(UserCtx *user, const char *fieldName)
Updates the local vector (including ghost points) from its corresponding global vector.
Definition setup.c:1361
UserCtx * user
Definition variables.h:528
PetscMPIInt rank
Definition variables.h:646
PetscInt block_number
Definition variables.h:712
SimCtx * simCtx
Back-pointer to the master simulation context.
Definition variables.h:814
Vec Centz
Definition variables.h:859
PetscReal Min_X
Definition variables.h:821
PetscInt KM
Definition variables.h:820
Vec lZet
Definition variables.h:858
UserMG usermg
Definition variables.h:764
PetscReal ren
Definition variables.h:691
PetscInt _this
Definition variables.h:824
PetscReal ry
Definition variables.h:825
PetscReal Max_Y
Definition variables.h:821
Vec Ucont
Definition variables.h:837
Vec Ubcs
Physical Cartesian velocity at boundary faces. Full 3D array but only boundary-face entries are meani...
Definition variables.h:121
PetscScalar x
Definition variables.h:101
Vec Centx
Definition variables.h:859
BCS Bcs
Definition variables.h:832
VerificationScalarConfig verificationScalar
Definition variables.h:700
PetscReal rz
Definition variables.h:825
Vec lCsi
Definition variables.h:858
Cmpnts AnalyticalUniformVelocity
Definition variables.h:698
PetscScalar z
Definition variables.h:101
Vec Ucat
Definition variables.h:837
PetscInt JM
Definition variables.h:820
PetscInt mglevels
Definition variables.h:535
PetscReal Min_Z
Definition variables.h:821
char AnalyticalSolutionType[PETSC_MAX_PATH_LEN]
Definition variables.h:676
PetscReal Max_X
Definition variables.h:821
PetscReal Min_Y
Definition variables.h:821
DMDALocalInfo info
Definition variables.h:818
PetscScalar y
Definition variables.h:101
PetscInt IM
Definition variables.h:820
Vec lEta
Definition variables.h:858
Vec Cent
Definition variables.h:858
MGCtx * mgctx
Definition variables.h:538
Vec Centy
Definition variables.h:859
PetscReal rx
Definition variables.h:825
PetscReal ti
Definition variables.h:652
PetscReal Max_Z
Definition variables.h:821
A 3D point or vector with PetscScalar components.
Definition variables.h:100
The master context for the entire simulation.
Definition variables.h:643
User-defined context containing data specific to a single computational grid level.
Definition variables.h:811
Verification-only analytical scalar override settings.
Definition variables.h:207