PICurv 0.1.0
A Parallel Particle-In-Cell Solver for Curvilinear LES
Loading...
Searching...
No Matches
interpolation.h
Go to the documentation of this file.
1#ifndef INTERPOLATION_H
2#define INTERPOLATION_H
3
4#include <petscpf.h>
5#include <petscdmswarm.h>
6#include <stdlib.h>
7#include <time.h>
8#include <math.h>
9#include <petsctime.h>
10#include <petscdmcomposite.h>
11#include <petscerror.h>
12#include <petscsys.h>
13
14// Include additional headers
15#include "variables.h" // Shared type definitions
16#include "ParticleSwarm.h" // Particle swarm functions
17#include "walkingsearch.h" // Particle location functions
18#include "grid.h" // Grid functions
19#include "logging.h" // Logging macros
20#include "io.h" // Data Input and Output functions
21#include "setup.h" // Setup functions that are used across the codebase
22
23// Macros and constants
24#define NUM_WEIGHTS 8 // Number of weights in trilinear interpolation
25
26/**
27 * @brief Generic macro to call the appropriate interpolation function based on the field type.
28 *
29 * This macro will select either the scalar or vector interpolation function based on
30 * the type of the 'field' pointer. It also performs a compile-time check that 'centfield'
31 * is of the same type as 'field'. If the types are not the same, a compile-time error is produced.
32 *
33 * Usage:
34 * InterpolateFieldFromCornerToCenter(field, centfield, user);
35 */
36#define InterpolateFieldFromCornerToCenter(field, centfield, user) \
37 ( (void)sizeof(char[1 - 2*!!(!__builtin_types_compatible_p(typeof(field), typeof(centfield)))]),\
38 _Generic((field), \
39 PetscReal ***: InterpolateFieldFromCornerToCenter_Scalar, \
40 Cmpnts ***: InterpolateFieldFromCornerToCenter_Vector \
41 )(field, centfield, user) )
42
43/**
44 * @brief Macro to dispatch to the correct scalar or vector center-to-corner function
45 * based on a runtime block size variable.
46 *
47 * This macro uses a ternary operator to inspect the runtime value of 'blockSize' and
48 * select the appropriate implementation. It expects the input and output pointers to be void*
49 * and handles casting them to the correct, strongly-typed pointers.
50 *
51 * @param blockSize The runtime block size (1 for scalar, 3 for vector).
52 * @param centfield_ptr The input void* pointer to the 3D cell-centered data array.
53 * @param corner_ptr The output void* pointer to the 3D corner data array.
54 * @param user_ctx The UserCtx structure.
55 */
56#define InterpolateFieldFromCenterToCorner(blockSize, centfield_ptr, corner_ptr, user_ctx) \
57 ( (blockSize) == 1 ? \
58 InterpolateFieldFromCenterToCorner_Scalar((PetscReal***)(centfield_ptr), (PetscReal***)(corner_ptr), (user_ctx)) : \
59 InterpolateFieldFromCenterToCorner_Vector((Cmpnts***)(centfield_ptr), (Cmpnts***)(corner_ptr), (user_ctx)) \
60 )
61
62
63/**
64 * @brief A type-generic macro that interpolates a field from corner nodes to all face centers.
65 *
66 * This macro uses C11's _Generic feature to dispatch to the appropriate underlying
67 * function based on the type of the input `corner_arr`.
68 *
69 * - If `corner_arr` is of type `PetscReal***`, it calls `InterpolateCornerToFaceCenter_Scalar`.
70 * - If `corner_arr` is of type `Cmpnts***`, it calls `InterpolateCornerToFaceCenter_Vector`.
71 *
72 * @param corner_arr Ghosted node-centered array (global indexing). The type of this
73 * argument determines which function is called.
74 * @param faceX_arr Local array for X-faces.
75 * @param faceY_arr Local array for Y-faces.
76 * @param faceZ_arr Local array for Z-faces.
77 * @param user_ctx User context containing DMDA 'fda'.
78 */
79#define InterpolateCornerToFaceCenter(corner_arr, faceX_arr, faceY_arr, faceZ_arr, user_ctx) \
80 _Generic((corner_arr), \
81 PetscReal***: InterpolateCornerToFaceCenter_Scalar, \
82 Cmpnts***: InterpolateCornerToFaceCenter_Vector \
83 )(corner_arr, faceX_arr, faceY_arr, faceZ_arr, user_ctx)
84
85/**
86 * @brief Macro that calls either the scalar or vector piecewise interpolation function
87 * based on the type of the `fieldPtr` parameter (3D array).
88 *
89 * Usage example:
90 *
91 * // For scalar:
92 * PetscReal ***fieldScal;
93 * PetscReal outVal;
94 * PieceWiseLinearInterpolation(fieldName, fieldScal, i, j, k, &outVal);
95 *
96 * // For vector:
97 * Cmpnts ***fieldVec;
98 * Cmpnts vec;
99 * PieceWiseLinearInterpolation(fieldName, fieldVec, i, j, k, &vec);
100 */
101#define PieceWiseLinearInterpolation(fieldName, fieldPtr, i, j, k, outPtr) \
102 _Generic((fieldPtr), \
103 PetscReal ***: PieceWiseLinearInterpolation_Scalar, \
104 Cmpnts ***: PieceWiseLinearInterpolation_Vector \
105 )(fieldName, fieldPtr, i, j, k, outPtr)
106
107/**
108 * @brief Returns the scalar value from the input cell index without blending.
109 *
110 * This is a first-order piecewise-constant lookup helper used in workflows that
111 * intentionally disable trilinear blending.
112 *
113 * @param[in] fieldName Field label used for diagnostics.
114 * @param[in] fieldScal Scalar field array indexed as [k][j][i].
115 * @param[in] iCell Cell i-index.
116 * @param[in] jCell Cell j-index.
117 * @param[in] kCell Cell k-index.
118 * @param[out] val Output scalar value.
119 * @return PetscErrorCode 0 on success.
120 */
122 const char *fieldName,
123 PetscReal ***fieldScal,
124 PetscInt iCell,
125 PetscInt jCell,
126 PetscInt kCell,
127 PetscReal *val);
128
129/**
130 * @brief Returns the vector value from the input cell index without blending.
131 *
132 * This is a first-order piecewise-constant lookup helper used in workflows that
133 * intentionally disable trilinear blending.
134 *
135 * @param[in] fieldName Field label used for diagnostics.
136 * @param[in] fieldVec Vector field array indexed as [k][j][i].
137 * @param[in] iCell Cell i-index.
138 * @param[in] jCell Cell j-index.
139 * @param[in] kCell Cell k-index.
140 * @param[out] vec Output vector value.
141 * @return PetscErrorCode 0 on success.
142 */
144 const char *fieldName,
145 Cmpnts ***fieldVec,
146 PetscInt iCell,
147 PetscInt jCell,
148 PetscInt kCell,
149 Cmpnts *vec);
150
151/**
152 * @brief Macro that calls either the scalar or vector trilinear interpolation function
153 * based on the type of the `fieldPtr` parameter (3D array).
154 *
155 * Usage example:
156 *
157 * PetscReal result;
158 * Cmpnts vec;
159 *
160 * // For scalars:
161 * TrilinearInterpolation(fieldName, fieldScal, i, j, k, a1, a2, a3, &result);
162 *
163 * // For vectors:
164 * TrilinearInterpolation(fieldName, fieldVec, i, j, k, a1, a2, a3, &vec);
165 */
166#define TrilinearInterpolation(fieldName, fieldPtr, i, j, k, a1, a2, a3, outPtr) \
167 _Generic((fieldPtr), \
168 PetscReal ***: TrilinearInterpolation_Scalar, \
169 Cmpnts ***: TrilinearInterpolation_Vector \
170 )(fieldName, fieldPtr, i, j, k, a1, a2, a3, outPtr)
171
172// Function declarations
173
174/**
175 * @brief Computes the trilinear interpolated scalar at a given point.
176 *
177 * each cell a PetscReal.
178 * This function uses the standard 8-corner trilinear formula via `ComputeTrilinearWeights()`.
179 * If a different scheme is desired, implement a new function with the same interface.
180 *
181 * @param fieldName A
182 * @param fieldScal 3D
183 * @param i Integral
184 * @param j Integral
185 * @param k Integral
186 * @param a1 Normalized
187 * @param a2 Normalized
188 * @param a3 Normalized
189 * @param val Pointer
190 * @return PetscErrorCode 0 on success.
191 */
192PetscErrorCode TrilinearInterpolation_Scalar(
193 const char *fieldName,
194 PetscReal ***fieldScal,
195 PetscInt i,
196 PetscInt j,
197 PetscInt k,
198 PetscReal a1,
199 PetscReal a2,
200 PetscReal a3,
201 PetscReal *val);
202
203/**
204 * @brief Computes the trilinear interpolated vector (e.g., velocity) at a given point.
205 *
206 * each cell of type Cmpnts.
207 * This function uses the standard 8-corner trilinear formula via `ComputeTrilinearWeights()`.
208 * If a different scheme is desired, implement a new function with the same interface.
209 *
210 * @param fieldName A
211 * @param fieldVec 3D
212 * @param i Integral
213 * @param j Integral
214 * @param k Integral
215 * @param a1 Normalized
216 * @param a2 Normalized
217 * @param a3 Normalized
218 * @param vec Pointer
219 * @return PetscErrorCode 0 on success.
220 */
221PetscErrorCode TrilinearInterpolation_Vector(
222 const char *fieldName,
223 Cmpnts ***fieldVec,
224 PetscInt i,
225 PetscInt j,
226 PetscInt k,
227 PetscReal a1,
228 PetscReal a2,
229 PetscReal a3,
230 Cmpnts *vec);
231
232/**
233 * @brief Interpolates a cell-centered field (scalar or vector) onto DMSwarm particles,
234 * using a robust, PETSc-idiomatic two-stage process.
235 *
236 * This function first converts the cell-centered input data to corner-node data,
237 * storing this intermediate result in a PETSc Vec to correctly handle the communication
238 * of ghost-point information across parallel ranks. It then performs a final trilinear
239 * interpolation from the ghosted corner data to each particle's location.
240 *
241 * Workflow:
242 * 1. Create temporary PETSc Vecs (`cornerGlobal`, `cornerLocal`) to manage the
243 * intermediate corner-node data, using the existing nodal DMDA (`user->fda`).
244 * 2. Call a dispatch macro that uses the runtime block size (`bs`) to select the
245 * correct underlying center-to-corner function, writing results into `cornerGlobal`.
246 * 3. Perform a ghost-point exchange (`DMGlobalToLocal`) to transfer the boundary
247 * data from `cornerGlobal` into the ghost regions of `cornerLocal`.
248 * 4. Loop over all local particles. For each particle:
249 * a. Convert its global cell index to a local index relative to the ghosted array.
250 * b. Check if the particle's interpolation stencil is fully contained within
251 * the owned+ghost region. If not, log a warning and set the result to zero.
252 * c. Perform the final trilinear interpolation using the ghosted `cornerLocal` data.
253 * 5. Restore all PETSc objects to prevent memory leaks.
254 *
255 * @param[in] user User context with DMDA, DMSwarm, etc.
256 * @param[in] fieldLocal_cellCentered Ghosted local Vec with cell-centered data (e.g., Ucat).
257 * @param[in] fieldName Human-readable field name for logging (e.g., "Ucat").
258 * @param[in] swarmOutFieldName Name of the DMSwarm field where interpolation results go.
259 *
260 * @return PetscErrorCode 0 on success.
261 */
262PetscErrorCode InterpolateEulerFieldToSwarm(
263 UserCtx *user,
264 Vec fieldLocal_cellCentered,
265 const char *fieldName,
266 const char *swarmOutFieldName);
267
268
269/**
270 * @brief Interpolates all relevant fields from the DMDA to the DMSwarm.
271 *
272 * Currently, it interpolates:
273 * - user->Ucat (vector field) into the DMSwarm field "swarmVelocity".
274 *
275 * To add more fields, duplicate the call to InterpolateOneFieldOverSwarm and provide:
276 * - The global Vec for that field (e.g. user->Tcat for temperature),
277 * - A human-readable field name (for logging),
278 * - A DMSwarm output field name (e.g. "swarmTemperature").
279 *
280 * @param[in,out] user Pointer to a UserCtx containing:
281 * - user->da (DM for the grid),
282 * - user->swarm (DMSwarm for particles),
283 * - user->Ucat (Vec for the vector field),
284 * - possibly more fields like user->Tcat, user->Pcat, etc.
285 *
286 * @return PetscErrorCode Returns 0 on success, non-zero on failure.
287 */
288PetscErrorCode InterpolateAllFieldsToSwarm(UserCtx *user);
289
290/**
291 * @brief Interpolates particle velocities using trilinear interpolation.
292 *
293 * @param[in] user Pointer to the user-defined context containing grid and swarm information.
294 * @return PetscErrorCode Returns 0 on success, non-zero on failure.
295 */
296
298
299/**
300 * @brief Safely interpolate a scalar field from corner nodes (from the coordinate DM)
301 * to cell centers (from the cell-centered DM) using the provided UserCtx.
302 *
303 * For each cell center in the physical region of the cell-centered DM (fda), this function
304 * averages the up to 8 surrounding scalar values from the coordinate DM (da). On boundaries,
305 * where fewer corners are available, a partial average is computed.
306 *
307 * The coordinate DM (da) is built on corners (IM+1 x JM+1 x KM+1) while the cell-centered DM (fda)
308 * covers the physical cells (IM x JM x KM). Index offsets are adjusted via DMDAGetLocalInfo.
309 *
310 * @param[in] field_arr 3D array of corner-based scalar data (from user->da).
311 * @param[out] centfield_arr 3D array for interpolated cell-center scalar data (for user->fda).
312 * @param[in] user User context containing:
313 * - da : DM for the coordinate (corner) data.
314 * - fda : DM for the cell-centered data.
315 *
316 * @return PetscErrorCode 0 on success.
317 */
319 PetscReal ***field_arr,
320 PetscReal ***centfield_arr,
321 UserCtx *user);
322
323/**
324 * @brief Safely interpolate a vector field from corner nodes (from the coordinate DM)
325 * to cell centers (from the cell-centered DM) using the provided UserCtx.
326 *
327 * For each cell center in the physical region of the cell-centered DM (fda), this function
328 * averages the 8 surrounding corner values from the coordinate DM (da). The coordinate DM
329 * (da) is built on corners (IM+1 x JM+1 x KM+1) while the cell-centered DM (fda) covers
330 * the physical cells (IM x JM x KM). Index offsets are adjusted using DMDAGetLocalInfo.
331 *
332 * @param[in] field_arr 3D array of corner-based vector data (from user->da).
333 * @param[out] centfield_arr 3D array for interpolated cell-center vector data (for user->fda).
334 * @param[in] user User context containing:
335 * - da : DM for the coordinate (corner) data.
336 * - fda : DM for the cell-centered data.
337 *
338 * @return PetscErrorCode 0 on success.
339 */
341 Cmpnts ***field_arr,
342 Cmpnts ***centfield_arr,
343 UserCtx *user);
344
345/**
346 * @brief Tests the InterpolateFieldFromCornerToCenter function by reproducing the Cent vector.
347 *
348 * This function serves as a unit test. It performs the following steps:
349 * 1. Takes the corner-centered nodal coordinates (from DMGetCoordinatesLocal) as input.
350 * 2. Uses the `InterpolateFieldFromCornerToCenter` macro to interpolate these coordinates to
351 * the cell centers, storing the result in a new temporary vector.
352 * 3. Compares this new vector with the `user->Cent` vector, which is assumed to have been
353 * computed by `ComputeCellCentersAndSpacing` and serves as the ground truth.
354 * 4. A 2-norm of the difference is computed. If it is below a small tolerance, the test passes.
355 *
356 * @note This function should be called immediately after `ComputeCellCentersAndSpacing` has
357 * been successfully executed.
358 *
359 * @param user The UserCtx for a specific grid level.
360 * @return PetscErrorCode 0 on success, or a PETSc error code on failure.
361 */
362PetscErrorCode TestCornerToCenterInterpolation(UserCtx *user);
363
364
365/**
366 * @brief Interpolates a vector field from cell centers to corner nodes.
367 *
368 * This version is adapted to write directly into a ghosted local array obtained from DMDAVecGetArray(),
369 * which allows using GLOBAL indices for writing to the OWNED portion of the array.
370 *
371 * @param[in] centfield_arr Input: 3D array (ghosted) of cell-centered data, accessed via GLOBAL indices.
372 * @param[out] corner_arr Output: 3D array (ghosted) where interpolated node values are stored,
373 * also accessed via GLOBAL indices for the owned part.
374 * @param[in] user User context containing DMDA information.
375 *
376 * @return PetscErrorCode 0 on success.
377 */
379 Cmpnts ***centfield_arr, /* Input: Ghosted local array from Vec (read) */
380 Cmpnts ***corner_arr, /* Output: Ghosted local array from Vec (write) */
381 UserCtx *user);
382
383/**
384 * @brief Interpolates a scalar field from cell centers to corner nodes.
385 *
386 * This version is adapted to write directly into a ghosted local array obtained from DMDAVecGetArray(),
387 * which allows using GLOBAL indices for writing to the OWNED portion of the array.
388 *
389 * @param[in] centfield_arr Input: 3D array (ghosted) of scalar data at cell centers, accessed via GLOBAL indices.
390 * @param[out] corner_arr Output: 3D array (ghosted) where interpolated node values are stored,
391 * also accessed via GLOBAL indices for the owned part.
392 * @param[in] user User context containing DMDA information.
393 *
394 * @return PetscErrorCode 0 on success.
395 */
397 PetscReal ***centfield_arr, /* Input: Ghosted local array from Vec (read) */
398 PetscReal ***corner_arr, /* Output: Ghosted local array from Vec (write) */
399 UserCtx *user);
400
401/**
402 * @brief Determines the target Eulerian DM and expected DOF for scattering a given particle field.
403 * @ingroup scatter_module
404 *
405 * Based on hardcoded rules mapping particle field names to user context DMs (da/fda).
406 * This function encapsulates the policy of where different fields should be scattered.
407 *
408 * @param[in] user Pointer to the UserCtx containing da and fda.
409 * @param[in] particleFieldName Name of the particle field (e.g., "P", "Ucat").
410 * @param[out] targetDM Pointer to store the determined target DM (da or fda).
411 * @param[out] expected_dof Pointer to store the expected DOF (1 or 3) for this field.
412 *
413 * @return PetscErrorCode Returns 0 on success, PETSC_ERR_ARG_UNKNOWN if field name is not recognized,
414 * or PETSC_ERR_ARG_NULL for NULL inputs.
415 */
416PetscErrorCode GetScatterTargetInfo(UserCtx *user, const char *particleFieldName,
417 DM *targetDM, PetscInt *expected_dof);
418
419/**
420 * @brief Retrieves the persistent local vector (e.g., lPsi, lUcat) for a given field name.
421 *
422 * @param user User context containing the persistent vectors.
423 * @param fieldName Name of the field ("Psi", "Ucat", etc.).
424 * @param localVec Output pointer to the vector.
425 * @return PetscErrorCode
426 */
427PetscErrorCode GetPersistentLocalVector(UserCtx *user, const char *fieldName, Vec *localVec);
428
429/**
430 * @brief Accumulates a particle field (scalar or vector) into a target grid sum vector.
431 * @ingroup scatter_module_internal
432 *
433 * This function iterates through local particles, identifies their cell using the
434 * "DMSwarm_CellID" field, and adds the particle's field value (`particleFieldName`)
435 * to the corresponding cell location in the `gridSumVec`. It handles both scalar
436 * (DOF=1) and vector (DOF=3) fields automatically based on the DOF of `gridSumDM`.
437 *
438 * IMPORTANT: The caller must ensure `gridSumVec` is zeroed before calling this
439 * function if a fresh sum calculation is desired.
440 *
441 * @param[in] swarm The DMSwarm containing particles.
442 * @param[in] particleFieldName Name of the field on the particles (must match DOF).
443 * @param[in] gridSumDM The DMDA associated with `gridSumVec`. Its DOF determines
444 * how many components are accumulated.
445 * @param[in,out] gridSumVec The Vec (associated with `gridSumDM`) to accumulate sums into.
446 *
447 * @return PetscErrorCode 0 on success. Errors if fields don't exist or DMs are incompatible.
448 */
449PetscErrorCode AccumulateParticleField(DM swarm, const char *particleFieldName,
450 DM gridSumDM, Vec gridSumVec);
451
452/**
453 * @brief Normalizes a grid vector of sums by a grid vector of counts to produce an average.
454 * @ingroup scatter_module_internal
455 *
456 * Calculates avgVec[i] = sumVec[i] / countVec[i] for each component of each
457 * OWNED cell where countVec[i] > 0. Sets avgVec[i] = 0 otherwise.
458 * Handles both scalar (DOF=1) and vector (DOF=3) data fields based on `dataDM`.
459 * Uses basic `VecGetArray`/`VecGetArrayRead` and manual index calculation.
460 *
461 * @param[in] countDM The DMDA associated with `countVec` (must have DOF=1).
462 * @param[in] countVec The Vec containing particle counts per cell (read-only).
463 * @param[in] dataDM The DMDA associated with `sumVec` and `avgVec` (must have DOF=1 or DOF=3).
464 * @param[in] sumVec The Vec containing the accumulated sums per cell (read-only).
465 * @param[in,out] avgVec The Vec where the calculated averages will be stored (overwritten).
466 *
467 * @return PetscErrorCode 0 on success.
468 */
469PetscErrorCode NormalizeGridVectorByCount(DM countDM, Vec countVec,
470 DM dataDM, Vec sumVec, Vec avgVec);
471
472/**
473 * @brief Scatters a particle field (scalar or vector) to the corresponding Eulerian field average.
474 * @ingroup scatter_module
475 *
476 * This is the main user-facing function. It determines the target Eulerian DM
477 * based on the `particleFieldName`, validates the provided `eulerFieldAverageVec`
478 * against the target DM, and then orchestrates the scatter operation by calling
479 * the internal helper function `ScatterParticleFieldToEulerField_Internal`.
480 * The final averaged result is stored IN-PLACE in `eulerFieldAverageVec`.
481 *
482 * @param[in] user Pointer to UserCtx containing da, fda, swarm, ParticleCount.
483 * @param[in] particleFieldName Name of the field in the DMSwarm (e.g., "P", "Ucat").
484 * @param[in,out] eulerFieldAverageVec Pre-created Vec associated with the correct target DM
485 * (implicitly da or fda). Result stored here.
486 *
487 * @return PetscErrorCode 0 on success. Errors on NULL input, unrecognized field name,
488 * or incompatible target vector.
489 */
490PetscErrorCode ScatterParticleFieldToEulerField(UserCtx *user,
491 const char *particleFieldName,
492 Vec eulerFieldAverageVec);
493
494/**
495 * @brief Scatters a predefined set of particle fields to their corresponding Eulerian fields.
496 * @ingroup scatter_module
497 *
498 * This convenience function calls the unified `ScatterParticleFieldToEulerField`
499 * for a standard set of fields ("P", potentially others). It assumes the target
500 * Eulerian Vec objects (e.g., `user->P`, `user->Ucat`) exist in the UserCtx structure
501 * and are correctly associated with their respective DMs (`user->da` or `user->fda`).
502 * It zeros the target Vecs before scattering.
503 *
504 * @param[in,out] user Pointer to the UserCtx structure containing all required DMs,
505 * Vecs (`ParticleCount`, target Eulerian fields like `P`, `Ucat`), and `swarm`.
506 *
507 * @return PetscErrorCode 0 on success. Errors if prerequisites (like ParticleCount)
508 * are missing or if underlying scatter calls fail.
509 */
511
512/**
513 * @brief Interpolates a scalar field from corner nodes to all face centers.
514 *
515 * This routine computes the average of the four corner-node values
516 * defining each face of a hexahedral cell:
517 * - X-faces (perpendicular to X): face between (i-1,i) in X-dir
518 * - Y-faces (perpendicular to Y): face between (j-1,j) in Y-dir
519 * - Z-faces (perpendicular to Z): face between (k-1,k) in Z-dir
520 *
521 * @param[in] corner_arr Ghosted node-centered array (global indexing) from user->fda.
522 * @param[out] faceX_arr Local array for X-faces sized [zm][ym][xm+1].
523 * @param[out] faceY_arr Local array for Y-faces sized [zm][ym+1][xm].
524 * @param[out] faceZ_arr Local array for Z-faces sized [zm+1][ym][xm].
525 * @param[in] user User context containing DMDA 'fda' and GetOwnedCellRange.
526 *
527 * @return PetscErrorCode 0 on success, non-zero on failure.
528 */
530 PetscReal ***corner_arr,
531 PetscReal ***faceX_arr,
532 PetscReal ***faceY_arr,
533 PetscReal ***faceZ_arr,
534 UserCtx *user);
535
536/**
537 * @brief Interpolates a vector field from corner nodes to all face centers.
538 *
539 * Identical to the scalar version, except it averages each component of the
540 * Cmpnts struct at the four corner-nodes per face.
541 *
542 * @param[in] corner_arr Ghosted 3-component array (global node indices).
543 * @param[out] faceX_arr Local array of Cmpnts for X-faces sized [zm][ym][xm+1].
544 * @param[out] faceY_arr Local array of Cmpnts for Y-faces sized [zm][ym+1][xm].
545 * @param[out] faceZ_arr Local array of Cmpnts for Z-faces sized [zm+1][ym][xm].
546 * @param[in] user User context containing DMDA 'fda'.
547 *
548 * @return PetscErrorCode 0 on success.
549 */
551 Cmpnts ***corner_arr,
552 Cmpnts ***faceX_arr,
553 Cmpnts ***faceY_arr,
554 Cmpnts ***faceZ_arr,
555 UserCtx *user);
556
557#endif // INTERPOLATION_H
Header file for Particle Swarm management functions.
Public interface for grid, solver, and metric setup routines.
PetscErrorCode NormalizeGridVectorByCount(DM countDM, Vec countVec, DM dataDM, Vec sumVec, Vec avgVec)
Normalizes a grid vector of sums by a grid vector of counts to produce an average.
PetscErrorCode GetPersistentLocalVector(UserCtx *user, const char *fieldName, Vec *localVec)
Retrieves the persistent local vector (e.g., lPsi, lUcat) for a given field name.
PetscErrorCode AccumulateParticleField(DM swarm, const char *particleFieldName, DM gridSumDM, Vec gridSumVec)
Accumulates a particle field (scalar or vector) into a target grid sum vector.
PetscErrorCode GetScatterTargetInfo(UserCtx *user, const char *particleFieldName, DM *targetDM, PetscInt *expected_dof)
Determines the target Eulerian DM and expected DOF for scattering a given particle field.
PetscErrorCode ScatterAllParticleFieldsToEulerFields(UserCtx *user)
Scatters a predefined set of particle fields to their corresponding Eulerian fields.
PetscErrorCode ScatterParticleFieldToEulerField(UserCtx *user, const char *particleFieldName, Vec eulerFieldAverageVec)
Scatters a particle field (scalar or vector) to the corresponding Eulerian field average.
PetscErrorCode InterpolateFieldFromCornerToCenter_Vector(Cmpnts ***field_arr, Cmpnts ***centfield_arr, UserCtx *user)
Safely interpolate a vector field from corner nodes (from the coordinate DM) to cell centers (from th...
PetscErrorCode InterpolateCornerToFaceCenter_Vector(Cmpnts ***corner_arr, Cmpnts ***faceX_arr, Cmpnts ***faceY_arr, Cmpnts ***faceZ_arr, UserCtx *user)
Interpolates a vector field from corner nodes to all face centers.
PetscErrorCode InterpolateFieldFromCornerToCenter_Scalar(PetscReal ***field_arr, PetscReal ***centfield_arr, UserCtx *user)
Safely interpolate a scalar field from corner nodes (from the coordinate DM) to cell centers (from th...
PetscErrorCode InterpolateAllFieldsToSwarm(UserCtx *user)
Interpolates all relevant fields from the DMDA to the DMSwarm.
PetscErrorCode TrilinearInterpolation_Vector(const char *fieldName, Cmpnts ***fieldVec, PetscInt i, PetscInt j, PetscInt k, PetscReal a1, PetscReal a2, PetscReal a3, Cmpnts *vec)
Computes the trilinear interpolated vector (e.g., velocity) at a given point.
PetscErrorCode InterpolateParticleVelocities(UserCtx *user)
Interpolates particle velocities using trilinear interpolation.
PetscErrorCode TrilinearInterpolation_Scalar(const char *fieldName, PetscReal ***fieldScal, PetscInt i, PetscInt j, PetscInt k, PetscReal a1, PetscReal a2, PetscReal a3, PetscReal *val)
Computes the trilinear interpolated scalar at a given point.
PetscErrorCode PieceWiseLinearInterpolation_Scalar(const char *fieldName, PetscReal ***fieldScal, PetscInt iCell, PetscInt jCell, PetscInt kCell, PetscReal *val)
Returns the scalar value from the input cell index without blending.
PetscErrorCode TestCornerToCenterInterpolation(UserCtx *user)
Tests the InterpolateFieldFromCornerToCenter function by reproducing the Cent vector.
PetscErrorCode InterpolateFieldFromCenterToCorner_Vector(Cmpnts ***centfield_arr, Cmpnts ***corner_arr, UserCtx *user)
Interpolates a vector field from cell centers to corner nodes.
PetscErrorCode InterpolateCornerToFaceCenter_Scalar(PetscReal ***corner_arr, PetscReal ***faceX_arr, PetscReal ***faceY_arr, PetscReal ***faceZ_arr, UserCtx *user)
Interpolates a scalar field from corner nodes to all face centers.
PetscErrorCode InterpolateEulerFieldToSwarm(UserCtx *user, Vec fieldLocal_cellCentered, const char *fieldName, const char *swarmOutFieldName)
Interpolates a cell-centered field (scalar or vector) onto DMSwarm particles, using a robust,...
PetscErrorCode PieceWiseLinearInterpolation_Vector(const char *fieldName, Cmpnts ***fieldVec, PetscInt iCell, PetscInt jCell, PetscInt kCell, Cmpnts *vec)
Returns the vector value from the input cell index without blending.
PetscErrorCode InterpolateFieldFromCenterToCorner_Scalar(PetscReal ***centfield_arr, PetscReal ***corner_arr, UserCtx *user)
Interpolates a scalar field from cell centers to corner nodes.
Public interface for data input/output routines.
Logging utilities and macros for PETSc-based applications.
Main header file for a complex fluid dynamics solver.
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:811
Header file for particle location functions using the walking search algorithm.