PICurv 0.1.0
A Parallel Particle-In-Cell Solver for Curvilinear LES
Loading...
Searching...
No Matches
Macros | Functions
interpolation.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 <petscerror.h>
#include <petscsys.h>
#include "variables.h"
#include "ParticleSwarm.h"
#include "walkingsearch.h"
#include "grid.h"
#include "logging.h"
#include "io.h"
#include "setup.h"
Include dependency graph for interpolation.h:
This graph shows which files directly or indirectly include this file:

Go to the source code of this file.

Macros

#define NUM_WEIGHTS   8
 
#define InterpolateFieldFromCornerToCenter(field, centfield, user)
 Generic macro to call the appropriate interpolation function based on the field type.
 
#define InterpolateFieldFromCenterToCorner(blockSize, centfield_ptr, corner_ptr, user_ctx)
 Macro to dispatch to the correct scalar or vector center-to-corner function based on a runtime block size variable.
 
#define InterpolateCornerToFaceCenter(corner_arr, faceX_arr, faceY_arr, faceZ_arr, user_ctx)
 A type-generic macro that interpolates a field from corner nodes to all face centers.
 
#define PieceWiseLinearInterpolation(fieldName, fieldPtr, i, j, k, outPtr)
 Macro that calls either the scalar or vector piecewise interpolation function based on the type of the fieldPtr parameter (3D array).
 
#define TrilinearInterpolation(fieldName, fieldPtr, i, j, k, a1, a2, a3, outPtr)
 Macro that calls either the scalar or vector trilinear interpolation function based on the type of the fieldPtr parameter (3D array).
 

Functions

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 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 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 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 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, PETSc-idiomatic two-stage process.
 
PetscErrorCode InterpolateAllFieldsToSwarm (UserCtx *user)
 Interpolates all relevant fields from the DMDA to the DMSwarm.
 
PetscErrorCode InterpolateParticleVelocities (UserCtx *user)
 Interpolates particle velocities using trilinear interpolation.
 
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 the cell-centered DM) using the provided UserCtx.
 
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 the cell-centered DM) using the provided UserCtx.
 
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 InterpolateFieldFromCenterToCorner_Scalar (PetscReal ***centfield_arr, PetscReal ***corner_arr, UserCtx *user)
 Interpolates a scalar field from cell centers to corner nodes.
 
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 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 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 ScatterParticleFieldToEulerField (UserCtx *user, const char *particleFieldName, Vec eulerFieldAverageVec)
 Scatters a particle field (scalar or vector) to the corresponding Eulerian field average.
 
PetscErrorCode ScatterAllParticleFieldsToEulerFields (UserCtx *user)
 Scatters a predefined set of particle fields to their corresponding Eulerian fields.
 
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 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.
 

Macro Definition Documentation

◆ NUM_WEIGHTS

#define NUM_WEIGHTS   8

Definition at line 24 of file interpolation.h.

◆ InterpolateFieldFromCornerToCenter

#define InterpolateFieldFromCornerToCenter (   field,
  centfield,
  user 
)
Value:
( (void)sizeof(char[1 - 2*!!(!__builtin_types_compatible_p(typeof(field), typeof(centfield)))]),\
_Generic((field), \
)(field, centfield, user) )
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 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...
A 3D point or vector with PetscScalar components.
Definition variables.h:100

Generic macro to call the appropriate interpolation function based on the field type.

This macro will select either the scalar or vector interpolation function based on the type of the 'field' pointer. It also performs a compile-time check that 'centfield' is of the same type as 'field'. If the types are not the same, a compile-time error is produced.

Usage: InterpolateFieldFromCornerToCenter(field, centfield, user);

Definition at line 36 of file interpolation.h.

◆ InterpolateFieldFromCenterToCorner

#define InterpolateFieldFromCenterToCorner (   blockSize,
  centfield_ptr,
  corner_ptr,
  user_ctx 
)
Value:
( (blockSize) == 1 ? \
InterpolateFieldFromCenterToCorner_Scalar((PetscReal***)(centfield_ptr), (PetscReal***)(corner_ptr), (user_ctx)) : \
InterpolateFieldFromCenterToCorner_Vector((Cmpnts***)(centfield_ptr), (Cmpnts***)(corner_ptr), (user_ctx)) \
)
PetscErrorCode InterpolateFieldFromCenterToCorner_Vector(Cmpnts ***centfield_arr, Cmpnts ***corner_arr, UserCtx *user)
Interpolates a vector field from cell centers to corner nodes.
PetscErrorCode InterpolateFieldFromCenterToCorner_Scalar(PetscReal ***centfield_arr, PetscReal ***corner_arr, UserCtx *user)
Interpolates a scalar field from cell centers to corner nodes.

Macro to dispatch to the correct scalar or vector center-to-corner function based on a runtime block size variable.

This macro uses a ternary operator to inspect the runtime value of 'blockSize' and select the appropriate implementation. It expects the input and output pointers to be void* and handles casting them to the correct, strongly-typed pointers.

Parameters
blockSizeThe runtime block size (1 for scalar, 3 for vector).
centfield_ptrThe input void* pointer to the 3D cell-centered data array.
corner_ptrThe output void* pointer to the 3D corner data array.
user_ctxThe UserCtx structure.

Definition at line 56 of file interpolation.h.

58 : \
59 InterpolateFieldFromCenterToCorner_Vector((Cmpnts***)(centfield_ptr), (Cmpnts***)(corner_ptr), (user_ctx)) \
60 )

◆ InterpolateCornerToFaceCenter

#define InterpolateCornerToFaceCenter (   corner_arr,
  faceX_arr,
  faceY_arr,
  faceZ_arr,
  user_ctx 
)
Value:
_Generic((corner_arr), \
)(corner_arr, faceX_arr, faceY_arr, faceZ_arr, user_ctx)
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 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.

A type-generic macro that interpolates a field from corner nodes to all face centers.

This macro uses C11's _Generic feature to dispatch to the appropriate underlying function based on the type of the input corner_arr.

  • If corner_arr is of type PetscReal***, it calls InterpolateCornerToFaceCenter_Scalar.
  • If corner_arr is of type Cmpnts***, it calls InterpolateCornerToFaceCenter_Vector.
Parameters
corner_arrGhosted node-centered array (global indexing). The type of this argument determines which function is called.
faceX_arrLocal array for X-faces.
faceY_arrLocal array for Y-faces.
faceZ_arrLocal array for Z-faces.
user_ctxUser context containing DMDA 'fda'.

Definition at line 79 of file interpolation.h.

83 )(corner_arr, faceX_arr, faceY_arr, faceZ_arr, user_ctx)

◆ PieceWiseLinearInterpolation

#define PieceWiseLinearInterpolation (   fieldName,
  fieldPtr,
  i,
  j,
  k,
  outPtr 
)
Value:
_Generic((fieldPtr), \
)(fieldName, fieldPtr, i, j, k, outPtr)
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 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.

Macro that calls either the scalar or vector piecewise interpolation function based on the type of the fieldPtr parameter (3D array).

Usage example:

// For scalar: PetscReal ***fieldScal; PetscReal outVal; PieceWiseLinearInterpolation(fieldName, fieldScal, i, j, k, &outVal);

// For vector: Cmpnts ***fieldVec; Cmpnts vec; PieceWiseLinearInterpolation(fieldName, fieldVec, i, j, k, &vec);

Definition at line 101 of file interpolation.h.

105 )(fieldName, fieldPtr, i, j, k, outPtr)

◆ TrilinearInterpolation

#define TrilinearInterpolation (   fieldName,
  fieldPtr,
  i,
  j,
  k,
  a1,
  a2,
  a3,
  outPtr 
)
Value:
_Generic((fieldPtr), \
)(fieldName, fieldPtr, i, j, k, a1, a2, a3, outPtr)
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 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.

Macro that calls either the scalar or vector trilinear interpolation function based on the type of the fieldPtr parameter (3D array).

Usage example:

PetscReal result; Cmpnts vec;

// For scalars: TrilinearInterpolation(fieldName, fieldScal, i, j, k, a1, a2, a3, &result);

// For vectors: TrilinearInterpolation(fieldName, fieldVec, i, j, k, a1, a2, a3, &vec);

Definition at line 166 of file interpolation.h.

170 )(fieldName, fieldPtr, i, j, k, a1, a2, a3, outPtr)

Function Documentation

◆ PieceWiseLinearInterpolation_Scalar()

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.

This is a first-order piecewise-constant lookup helper used in workflows that intentionally disable trilinear blending.

Parameters
[in]fieldNameField label used for diagnostics.
[in]fieldScalScalar field array indexed as [k][j][i].
[in]iCellCell i-index.
[in]jCellCell j-index.
[in]kCellCell k-index.
[out]valOutput scalar value.
Returns
PetscErrorCode 0 on success.

Returns the scalar value from the input cell index without blending.

Local to this translation unit.

Definition at line 462 of file interpolation.c.

469{
470 PetscFunctionBegin;
471 *val = fieldScal[kCell][jCell][iCell];
472
473 // Optional logging
475 "Field '%s' at (i=%d, j=%d, k=%d) => val=%.6f\n",
476 fieldName, iCell, jCell, kCell, *val);
477
478 PetscFunctionReturn(0);
479}
#define LOCAL
Logging scope definitions for controlling message output.
Definition logging.h:44
#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
@ LOG_VERBOSE
Extremely detailed logs, typically for development use only.
Definition logging.h:33

◆ PieceWiseLinearInterpolation_Vector()

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.

This is a first-order piecewise-constant lookup helper used in workflows that intentionally disable trilinear blending.

Parameters
[in]fieldNameField label used for diagnostics.
[in]fieldVecVector field array indexed as [k][j][i].
[in]iCellCell i-index.
[in]jCellCell j-index.
[in]kCellCell k-index.
[out]vecOutput vector value.
Returns
PetscErrorCode 0 on success.

Returns the vector value from the input cell index without blending.

Local to this translation unit.

Definition at line 487 of file interpolation.c.

494{
495 PetscFunctionBegin;
496 vec->x = fieldVec[kCell][jCell][iCell].x;
497 vec->y = fieldVec[kCell][jCell][iCell].y;
498 vec->z = fieldVec[kCell][jCell][iCell].z;
499
500 // Optional logging
502 "Field '%s' at (i=%d, j=%d, k=%d) => (x=%.6f, y=%.6f, z=%.6f)\n",
503 fieldName, iCell, jCell, kCell, vec->x, vec->y, vec->z);
504
505 PetscFunctionReturn(0);
506}
PetscScalar x
Definition variables.h:101
PetscScalar z
Definition variables.h:101
PetscScalar y
Definition variables.h:101

◆ TrilinearInterpolation_Scalar()

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.

each cell a PetscReal. This function uses the standard 8-corner trilinear formula via ComputeTrilinearWeights(). If a different scheme is desired, implement a new function with the same interface.

Parameters
fieldNameA
fieldScal3D
iIntegral
jIntegral
kIntegral
a1Normalized
a2Normalized
a3Normalized
valPointer
Returns
PetscErrorCode 0 on success.

Computes the trilinear interpolated scalar at a given point.

Local to this translation unit.

Definition at line 577 of file interpolation.c.

587{
588 PetscFunctionBegin; // PETSc macro for error/stack tracing
589
590 // Compute the 8 corner weights
591 PetscReal wcorner[8];
592 ComputeTrilinearWeights(a1, a2, a3, wcorner);
593
594 // Offsets for cell corners
595 PetscInt i1 = i + 1;
596 PetscInt j1 = j + 1;
597 PetscInt k1 = k + 1;
598
599 // Initialize the output scalar
600 PetscReal sum = 0.0;
601
602 // Corner 0 => (i, j, k)
603 sum += wcorner[0] * fieldScal[k ][j ][i ];
604 // Corner 1 => (i+1, j, k)
605 sum += wcorner[1] * fieldScal[k ][j ][i1];
606 // Corner 2 => (i, j+1, k)
607 sum += wcorner[2] * fieldScal[k ][j1][i ];
608 // Corner 3 => (i+1, j+1, k)
609 sum += wcorner[3] * fieldScal[k ][j1][i1];
610 // Corner 4 => (i, j, k+1)
611 sum += wcorner[4] * fieldScal[k1][j ][i ];
612 // Corner 5 => (i+1, j, k+1)
613 sum += wcorner[5] * fieldScal[k1][j ][i1];
614 // Corner 6 => (i, j+1, k+1)
615 sum += wcorner[6] * fieldScal[k1][j1][i ];
616 // Corner 7 => (i+1, j+1, k+1)
617 sum += wcorner[7] * fieldScal[k1][j1][i1];
618
619 *val = sum;
620
621 // Logging (optional)
623 "Field '%s' at (i=%d, j=%d, k=%d), "
624 "a1=%.6f, a2=%.6f, a3=%.6f -> val=%.6f.\n",
625 fieldName, i, j, k, a1, a2, a3, *val);
626
627 // LOG_ALLOW_SYNC(GLOBAL, LOG_INFO,
628 // "TrilinearInterpolation_Scalar: Completed interpolation for field '%s' across local cells.\n",
629 // fieldName);
630
631 PetscFunctionReturn(0);
632}
static void ComputeTrilinearWeights(PetscReal a1, PetscReal a2, PetscReal a3, PetscReal *w)
Internal helper implementation: ComputeTrilinearWeights().
Here is the call graph for this function:
Here is the caller graph for this function:

◆ TrilinearInterpolation_Vector()

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.

each cell of type Cmpnts. This function uses the standard 8-corner trilinear formula via ComputeTrilinearWeights(). If a different scheme is desired, implement a new function with the same interface.

Parameters
fieldNameA
fieldVec3D
iIntegral
jIntegral
kIntegral
a1Normalized
a2Normalized
a3Normalized
vecPointer
Returns
PetscErrorCode 0 on success.

Computes the trilinear interpolated vector (e.g., velocity) at a given point.

Local to this translation unit.

Definition at line 641 of file interpolation.c.

651{
652 PetscFunctionBegin; // PETSc macro for error/stack tracing
653
654 // Compute the 8 corner weights
655 PetscErrorCode ierr;
656 PetscReal wcorner[8];
657 PetscMPIInt rank;
658
659 ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
660
661 LOG_ALLOW(LOCAL,LOG_VERBOSE,"[Rank %d] Computing Trilinear Weights.\n",rank);
662 ComputeTrilinearWeights(a1, a2, a3, wcorner);
663
664 LOG_ALLOW(LOCAL,LOG_VERBOSE,"[Rank %d] Trilinear weights computed for local cell %d,%d,%d.\n",rank,i,j,k);
665
666 // For partial interpolation, we'll keep track of how many corners are valid
667 // and how much sum of weights is used. Then we do a final normalization.
668 PetscReal sumW = 0.0;
669 Cmpnts accum = {0.0, 0.0, 0.0};
670
671 // The eight corner indices, with their weights:
672 // corners: (i,j,k), (i+1,j,k), (i,j+1,k), (i+1,j+1,k), etc.
673 // We store them in an array to iterate cleanly.
674 const PetscInt cornerOffsets[8][3] = {
675 {0, 0, 0},
676 {1, 0, 0},
677 {0, 1, 0},
678 {1, 1, 0},
679 {0, 0, 1},
680 {1, 0, 1},
681 {0, 1, 1},
682 {1, 1, 1}
683 };
684
685 // Weighted partial sum
686 for (PetscInt c = 0; c < 8; c++) {
687 const PetscInt di = cornerOffsets[c][0];
688 const PetscInt dj = cornerOffsets[c][1];
689 const PetscInt dk = cornerOffsets[c][2];
690 PetscInt iC = i + di;
691 PetscInt jC = j + dj;
692 PetscInt kC = k + dk;
693
694 /*
695 // skip if out of domain
696 // (Assuming you know global domain is [0..mx), [0..my), [0..mz).)
697 if (iC < 0 || iC >= (PetscInt)userGlobalMx ||
698 jC < 0 || jC >= (PetscInt)userGlobalMy ||
699 kC < 0 || kC >= (PetscInt)userGlobalMz)
700 {
701 // skip this corner
702 continue;
703 }
704
705 */
706
707 LOG_ALLOW(LOCAL,LOG_VERBOSE,"[Rank %d] %s[%d][%d][%d] = (%.4f,%.4f,%.4f).\n",rank,fieldName,kC,jC,iC,fieldVec[kC][jC][iC].x,fieldVec[kC][jC][iC].y,fieldVec[kC][jC][iC].z);
708
709 // Otherwise, accumulate
710 accum.x += wcorner[c] * fieldVec[kC][jC][iC].x;
711 accum.y += wcorner[c] * fieldVec[kC][jC][iC].y;
712 accum.z += wcorner[c] * fieldVec[kC][jC][iC].z;
713 sumW += wcorner[c];
714 }
715
716 // If sumW=0 => out-of-range or zero weighting => set (0,0,0)
717 if (sumW > 1.0e-14) {
718 vec->x = accum.x / sumW;
719 vec->y = accum.y / sumW;
720 vec->z = accum.z / sumW;
721 } else {
722 vec->x = 0.0; vec->y = 0.0; vec->z = 0.0;
723 }
724
725 PetscFunctionReturn(0);
726}
Here is the call graph for this function:
Here is the caller graph for this function:

◆ InterpolateEulerFieldToSwarm()

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, PETSc-idiomatic two-stage process.

This function first converts the cell-centered input data to corner-node data, storing this intermediate result in a PETSc Vec to correctly handle the communication of ghost-point information across parallel ranks. It then performs a final trilinear interpolation from the ghosted corner data to each particle's location.

Workflow:

  1. Create temporary PETSc Vecs (cornerGlobal, cornerLocal) to manage the intermediate corner-node data, using the existing nodal DMDA (user->fda).
  2. Call a dispatch macro that uses the runtime block size (bs) to select the correct underlying center-to-corner function, writing results into cornerGlobal.
  3. Perform a ghost-point exchange (DMGlobalToLocal) to transfer the boundary data from cornerGlobal into the ghost regions of cornerLocal.
  4. Loop over all local particles. For each particle: a. Convert its global cell index to a local index relative to the ghosted array. b. Check if the particle's interpolation stencil is fully contained within the owned+ghost region. If not, log a warning and set the result to zero. c. Perform the final trilinear interpolation using the ghosted cornerLocal data.
  5. Restore all PETSc objects to prevent memory leaks.
Parameters
[in]userUser context with DMDA, DMSwarm, etc.
[in]fieldLocal_cellCenteredGhosted local Vec with cell-centered data (e.g., Ucat).
[in]fieldNameHuman-readable field name for logging (e.g., "Ucat").
[in]swarmOutFieldNameName of the DMSwarm field where interpolation results go.
Returns
PetscErrorCode 0 on success.

Interpolates a cell-centered field (scalar or vector) onto DMSwarm particles, using a robust, PETSc-idiomatic two-stage process.

Routes to InterpolateEulerFieldFromCenterToSwarm (direct trilinear, second-order) or InterpolateEulerFieldFromCornerToSwarm (corner-averaged, legacy) based on user->simCtx->interpolationMethod.

Definition at line 1371 of file interpolation.c.

1376{
1377 PetscErrorCode ierr;
1378 PetscFunctionBegin;
1379
1381 ierr = InterpolateEulerFieldFromCenterToSwarm(user, fieldLocal_cellCentered,
1382 fieldName, swarmOutFieldName); CHKERRQ(ierr);
1383 } else {
1384 ierr = InterpolateEulerFieldFromCornerToSwarm(user, fieldLocal_cellCentered,
1385 fieldName, swarmOutFieldName); CHKERRQ(ierr);
1386 }
1387
1388 PetscFunctionReturn(0);
1389}
static PetscErrorCode InterpolateEulerFieldFromCornerToSwarm(UserCtx *user, Vec fieldLocal_cellCentered, const char *fieldName, const char *swarmOutFieldName)
Corner-averaged interpolation path (legacy).
static PetscErrorCode InterpolateEulerFieldFromCenterToSwarm(UserCtx *user, Vec fieldLocal_cellCentered, const char *fieldName, const char *swarmOutFieldName)
Direct cell-center trilinear interpolation (second-order on curvilinear grids).
SimCtx * simCtx
Back-pointer to the master simulation context.
Definition variables.h:814
InterpolationMethod interpolationMethod
Definition variables.h:744
@ INTERP_TRILINEAR
Definition variables.h:522
Here is the call graph for this function:
Here is the caller graph for this function:

◆ InterpolateAllFieldsToSwarm()

PetscErrorCode InterpolateAllFieldsToSwarm ( UserCtx user)

Interpolates all relevant fields from the DMDA to the DMSwarm.

Currently, it interpolates:

  • user->Ucat (vector field) into the DMSwarm field "swarmVelocity".

To add more fields, duplicate the call to InterpolateOneFieldOverSwarm and provide:

  • The global Vec for that field (e.g. user->Tcat for temperature),
  • A human-readable field name (for logging),
  • A DMSwarm output field name (e.g. "swarmTemperature").
Parameters
[in,out]userPointer to a UserCtx containing:
  • user->da (DM for the grid),
  • user->swarm (DMSwarm for particles),
  • user->Ucat (Vec for the vector field),
  • possibly more fields like user->Tcat, user->Pcat, etc.
Returns
PetscErrorCode Returns 0 on success, non-zero on failure.

Interpolates all relevant fields from the DMDA to the DMSwarm.

Local to this translation unit.

Definition at line 1397 of file interpolation.c.

1398{
1399 PetscErrorCode ierr;
1400 PetscMPIInt rank;
1401 PetscFunctionBegin;
1402
1404
1405 ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank); CHKERRQ(ierr);
1406
1407 /*
1408 1) Interpolate the 'velocity' field (user->Ucat) into DMSwarm's "swarmVelocity".
1409 - The function InterpolateOneFieldOverSwarm() loops over all local particles,
1410 retrieves iCell, jCell, kCell, (a1,a2,a3) from the DMSwarm,
1411 and calls the appropriate trilinear interpolation routine.
1412
1413 - "velocity" => just a label used in logging or debugging.
1414 - "swarmVelocity" => the DMSwarm field name where we store the result
1415 (assume dof=3 if user->Ucat has blockSize=3).
1416 */
1417
1419 " Interpolation of ucat to velocity begins on rank %d.\n",rank);
1420 // Make sure to pass the *LOCAL* Vector to the function below!
1421 ierr = InterpolateEulerFieldToSwarm(user, user->lUcat,
1422 "Ucat",
1423 "velocity"); CHKERRQ(ierr);
1425 "Diffusivity", "Diffusivity"); CHKERRQ(ierr);
1426
1428 "DiffusivityGradient", "DiffusivityGradient"); CHKERRQ(ierr);
1429 /* Add fields as necessary here*/
1430
1431 ierr = MPI_Barrier(PETSC_COMM_WORLD); CHKERRQ(ierr);
1433 "[rank %d]Completed Interpolateting all fields to the swarm.\n",rank);
1434
1436
1437 PetscFunctionReturn(0);
1438}
PetscErrorCode InterpolateEulerFieldToSwarm(UserCtx *user, Vec fieldLocal_cellCentered, const char *fieldName, const char *swarmOutFieldName)
Dispatches grid-to-particle interpolation to the method selected in the control file.
#define PROFILE_FUNCTION_END
Marks the end of a profiled code block.
Definition logging.h:779
@ 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
Vec lDiffusivityGradient
Definition variables.h:841
Vec lUcat
Definition variables.h:837
Vec lDiffusivity
Definition variables.h:840
Here is the call graph for this function:
Here is the caller graph for this function:

◆ InterpolateParticleVelocities()

PetscErrorCode InterpolateParticleVelocities ( UserCtx user)

Interpolates particle velocities using trilinear interpolation.

Parameters
[in]userPointer to the user-defined context containing grid and swarm information.
Returns
PetscErrorCode Returns 0 on success, non-zero on failure.

◆ InterpolateFieldFromCornerToCenter_Scalar()

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 the cell-centered DM) using the provided UserCtx.

For each cell center in the physical region of the cell-centered DM (fda), this function averages the up to 8 surrounding scalar values from the coordinate DM (da). On boundaries, where fewer corners are available, a partial average is computed.

The coordinate DM (da) is built on corners (IM+1 x JM+1 x KM+1) while the cell-centered DM (fda) covers the physical cells (IM x JM x KM). Index offsets are adjusted via DMDAGetLocalInfo.

Parameters
[in]field_arr3D array of corner-based scalar data (from user->da).
[out]centfield_arr3D array for interpolated cell-center scalar data (for user->fda).
[in]userUser context containing:
  • da : DM for the coordinate (corner) data.
  • fda : DM for the cell-centered data.
Returns
PetscErrorCode 0 on success.

Safely interpolate a scalar field from corner nodes (from the coordinate DM) to cell centers (from the cell-centered DM) using the provided UserCtx.

Local to this translation unit.

Definition at line 86 of file interpolation.c.

90{
91 PetscErrorCode ierr;
92 DMDALocalInfo info;
93
94 PetscFunctionBeginUser;
95
96 ierr = DMDAGetLocalInfo(user->da, &info); CHKERRQ(ierr);
97
98 // Get local and global grid dimensions
99 PetscInt xs = info.xs, xe = info.xs + info.xm;
100 PetscInt ys = info.ys, ye = info.ys + info.ym;
101 PetscInt zs = info.zs, ze = info.zs + info.zm;
102 PetscInt mx = info.mx, my = info.my, mz = info.mz;
103
104 // Determine loop bounds to compute for interior cells only, matching the code's convention.
105 // Start at index 1 if the process owns the global boundary at index 0.
106 PetscInt is = (xs == 0) ? 1 : xs;
107 PetscInt js = (ys == 0) ? 1 : ys;
108 PetscInt ks = (zs == 0) ? 1 : zs;
109
110 // Stop one cell short if the process owns the global boundary at the max index.
111 PetscInt ie = (xe == mx) ? xe - 1 : xe;
112 PetscInt je = (ye == my) ? ye - 1 : ye;
113 PetscInt ke = (ze == mz) ? ze - 1 : ze;
114
115 // Loop over the locally owned INTERIOR cells.
116 for (PetscInt k = ks; k < ke; k++) {
117 for (PetscInt j = js; j < je; j++) {
118 for (PetscInt i = is; i < ie; i++) {
119 // Calculate cell center value as the average of its 8 corner nodes
120 centfield_arr[k][j][i] = 0.125 * (field_arr[k][j][i] + field_arr[k][j-1][i] +
121 field_arr[k-1][j][i] + field_arr[k-1][j-1][i] +
122 field_arr[k][j][i-1] + field_arr[k][j-1][i-1] +
123 field_arr[k-1][j][i-1] + field_arr[k-1][j-1][i-1]);
124 }
125 }
126 }
127
128 PetscFunctionReturn(0);
129}
Here is the caller graph for this function:

◆ InterpolateFieldFromCornerToCenter_Vector()

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 the cell-centered DM) using the provided UserCtx.

For each cell center in the physical region of the cell-centered DM (fda), this function averages the 8 surrounding corner values from the coordinate DM (da). The coordinate DM (da) is built on corners (IM+1 x JM+1 x KM+1) while the cell-centered DM (fda) covers the physical cells (IM x JM x KM). Index offsets are adjusted using DMDAGetLocalInfo.

Parameters
[in]field_arr3D array of corner-based vector data (from user->da).
[out]centfield_arr3D array for interpolated cell-center vector data (for user->fda).
[in]userUser context containing:
  • da : DM for the coordinate (corner) data.
  • fda : DM for the cell-centered data.
Returns
PetscErrorCode 0 on success.

Safely interpolate a vector field from corner nodes (from the coordinate DM) to cell centers (from the cell-centered DM) using the provided UserCtx.

Local to this translation unit.

Definition at line 25 of file interpolation.c.

29{
30 PetscErrorCode ierr;
31 DMDALocalInfo info;
32
33 PetscFunctionBeginUser;
34
35 ierr = DMDAGetLocalInfo(user->fda, &info); CHKERRQ(ierr);
36
37 // Get local and global grid dimensions
38 PetscInt xs = info.xs, xe = info.xs + info.xm;
39 PetscInt ys = info.ys, ye = info.ys + info.ym;
40 PetscInt zs = info.zs, ze = info.zs + info.zm;
41 PetscInt mx = info.mx, my = info.my, mz = info.mz;
42
43 // Determine loop bounds to compute for interior cells only, matching the code's convention.
44 // Start at index 1 if the process owns the global boundary at index 0.
45 PetscInt is = (xs == 0) ? 1 : xs;
46 PetscInt js = (ys == 0) ? 1 : ys;
47 PetscInt ks = (zs == 0) ? 1 : zs;
48
49 // Stop one cell short if the process owns the global boundary at the max index.
50 PetscInt ie = (xe == mx) ? xe - 1 : xe;
51 PetscInt je = (ye == my) ? ye - 1 : ye;
52 PetscInt ke = (ze == mz) ? ze - 1 : ze;
53
54 // Loop over the locally owned INTERIOR cells.
55 for (PetscInt k = ks; k < ke; k++) {
56 for (PetscInt j = js; j < je; j++) {
57 for (PetscInt i = is; i < ie; i++) {
58 // Calculate cell center value as the average of its 8 corner nodes
59 centfield_arr[k][j][i].x = 0.125 * (field_arr[k][j][i].x + field_arr[k][j-1][i].x +
60 field_arr[k-1][j][i].x + field_arr[k-1][j-1][i].x +
61 field_arr[k][j][i-1].x + field_arr[k][j-1][i-1].x +
62 field_arr[k-1][j][i-1].x + field_arr[k-1][j-1][i-1].x);
63
64 centfield_arr[k][j][i].y = 0.125 * (field_arr[k][j][i].y + field_arr[k][j-1][i].y +
65 field_arr[k-1][j][i].y + field_arr[k-1][j-1][i].y +
66 field_arr[k][j][i-1].y + field_arr[k][j-1][i-1].y +
67 field_arr[k-1][j][i-1].y + field_arr[k-1][j-1][i-1].y);
68
69 centfield_arr[k][j][i].z = 0.125 * (field_arr[k][j][i].z + field_arr[k][j-1][i].z +
70 field_arr[k-1][j][i].z + field_arr[k-1][j-1][i].z +
71 field_arr[k][j][i-1].z + field_arr[k][j-1][i-1].z +
72 field_arr[k-1][j][i-1].z + field_arr[k-1][j-1][i-1].z);
73 }
74 }
75 }
76
77 PetscFunctionReturn(0);
78}
Here is the caller graph for this function:

◆ TestCornerToCenterInterpolation()

PetscErrorCode TestCornerToCenterInterpolation ( UserCtx user)

Tests the InterpolateFieldFromCornerToCenter function by reproducing the Cent vector.

This function serves as a unit test. It performs the following steps:

  1. Takes the corner-centered nodal coordinates (from DMGetCoordinatesLocal) as input.
  2. Uses the InterpolateFieldFromCornerToCenter macro to interpolate these coordinates to the cell centers, storing the result in a new temporary vector.
  3. Compares this new vector with the user->Cent vector, which is assumed to have been computed by ComputeCellCentersAndSpacing and serves as the ground truth.
  4. A 2-norm of the difference is computed. If it is below a small tolerance, the test passes.
Note
This function should be called immediately after ComputeCellCentersAndSpacing has been successfully executed.
Parameters
userThe UserCtx for a specific grid level.
Returns
PetscErrorCode 0 on success, or a PETSc error code on failure.

Tests the InterpolateFieldFromCornerToCenter function by reproducing the Cent vector.

Local to this translation unit.

Definition at line 137 of file interpolation.c.

138{
139 PetscErrorCode ierr;
140 Vec lCoords, TestCent;
141 Cmpnts ***coor_arr, ***test_cent_arr;
142 PetscReal diff_norm;
143
144 PetscFunctionBeginUser;
145
146 // 1. Create a temporary vector to hold the result of our interpolation.
147 // It must have the same layout and size as the ground-truth user->Cent vector.
148 ierr = VecDuplicate(user->Cent, &TestCent); CHKERRQ(ierr);
149
150 // 2. Get the input (corner coordinates) and output (our test vector) arrays.
151 ierr = DMGetCoordinatesLocal(user->da, &lCoords); CHKERRQ(ierr);
152 ierr = DMDAVecGetArrayRead(user->fda, lCoords, &coor_arr); CHKERRQ(ierr);
153 ierr = DMDAVecGetArray(user->fda, TestCent, &test_cent_arr); CHKERRQ(ierr);
154
155 // 3. Call the generic interpolation macro.
156 // The macro will see that `coor_arr` is of type `Cmpnts***` and correctly
157 // call the `InterpolateFieldFromCornerToCenter_Vector` function.
158 ierr = InterpolateFieldFromCornerToCenter(coor_arr, test_cent_arr, user); CHKERRQ(ierr);
159
160 // 4. Restore the arrays.
161 ierr = DMDAVecRestoreArrayRead(user->fda, lCoords, &coor_arr); CHKERRQ(ierr);
162 ierr = DMDAVecRestoreArray(user->fda, TestCent, &test_cent_arr); CHKERRQ(ierr);
163
164 // 5. IMPORTANT: Assemble the vector so its values are communicated across processors
165 // and it's ready for global operations like VecNorm.
166 ierr = VecAssemblyBegin(TestCent); CHKERRQ(ierr);
167 ierr = VecAssemblyEnd(TestCent); CHKERRQ(ierr);
168
169 // 6. Compare the result with the ground truth.
170 // We compute TestCent = -1.0 * user->Cent + 1.0 * TestCent.
171 // This calculates the difference vector: TestCent - user->Cent.
172 ierr = VecAXPY(TestCent, -1.0, user->Cent); CHKERRQ(ierr);
173
174 // Now, compute the L2 norm of the difference vector. If the functions are
175 // identical, the norm should be zero (or very close due to floating point).
176 ierr = VecNorm(TestCent, NORM_2, &diff_norm); CHKERRQ(ierr);
177
178 // 7. Report the result and clean up.
179 if (diff_norm < 1.0e-12) {
180 LOG_ALLOW(GLOBAL,LOG_DEBUG,"[SUCCESS] Test passed. Norm of difference is %g.\n", (double)diff_norm);
181 } else {
182 LOG_ALLOW(GLOBAL,LOG_DEBUG, "[FAILURE] Test failed. Norm of difference is %g.\n", (double)diff_norm);
183 }
184
185 ierr = VecDestroy(&TestCent); CHKERRQ(ierr);
186
187 PetscFunctionReturn(0);
188}
#define InterpolateFieldFromCornerToCenter(field, centfield, user)
Generic macro to call the appropriate interpolation function based on the field type.
#define GLOBAL
Scope for global logging across all processes.
Definition logging.h:45
Vec Cent
Definition variables.h:858

◆ InterpolateFieldFromCenterToCorner_Vector()

PetscErrorCode InterpolateFieldFromCenterToCorner_Vector ( Cmpnts ***  centfield_arr,
Cmpnts ***  corner_arr,
UserCtx user 
)

Interpolates a vector field from cell centers to corner nodes.

This version is adapted to write directly into a ghosted local array obtained from DMDAVecGetArray(), which allows using GLOBAL indices for writing to the OWNED portion of the array.

Parameters
[in]centfield_arrInput: 3D array (ghosted) of cell-centered data, accessed via GLOBAL indices.
[out]corner_arrOutput: 3D array (ghosted) where interpolated node values are stored, also accessed via GLOBAL indices for the owned part.
[in]userUser context containing DMDA information.
Returns
PetscErrorCode 0 on success.

Interpolates a vector field from cell centers to corner nodes.

Local to this translation unit.

Definition at line 196 of file interpolation.c.

200{
201 PetscErrorCode ierr;
202 DMDALocalInfo info;
203 PetscMPIInt rank;
205 ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
206 ierr = DMDAGetLocalInfo(user->fda, &info); CHKERRQ(ierr);
207
208 // Node ownership range (GLOBAL indices)
209 PetscInt xs_node = info.xs, xm_node = info.xm, xe_node = xs_node + xm_node;
210 PetscInt ys_node = info.ys, ym_node = info.ym, ye_node = ys_node + ym_node;
211 PetscInt zs_node = info.zs, zm_node = info.zm, ze_node = zs_node + zm_node;
212
213 PetscInt nCellsX = info.mx - 2; // Number of cells in x-direction
214 PetscInt nCellsY = info.my - 2; // Number of cells in y-direction
215 PetscInt nCellsZ = info.mz - 2; // Number of cells in z-direction
216
217
218 // Global grid dimensions (used for valid cell check)
219 PetscInt IM = info.mx - 1; // Total nodes in i-direction
220 PetscInt JM = info.my - 1; // Total nodes in j-direction
221 PetscInt KM = info.mz - 1; // Total nodes in k-direction
222
224 "[Rank %d] Starting -- Node ownership k=%d..%d, j=%d..%d, i=%d..%d\n",
225 rank, zs_node, ze_node-1, ys_node, ye_node-1, xs_node, xe_node-1);
226
227 // Loop over the GLOBAL indices of the NODES owned by this processor
228 for (PetscInt k = zs_node; k < ze_node; k++) {
229 for (PetscInt j = ys_node; j < ye_node; j++) {
230 for (PetscInt i = xs_node; i < xe_node; i++) {
231 Cmpnts sum = {0.0, 0.0, 0.0};
232 PetscInt count = 0;
233
234 // DEBUG 1 TEST
235 /*
236 if(rank == 1 && i == 24 && j == 12 && k == 49){
237 PetscInt i_cell_A = i - 1;
238 PetscInt i_cell_B = i;
239
240 Cmpnts ucat_A = centfield_arr[k][j][i_cell_A]; // 23
241 Cmpnts ucat_B = centfield_arr[k][j][i_cell_B]; // 24 (out-of-bounds if IM=25)
242
243 PetscPrintf(PETSC_COMM_WORLD,"[Rank %d] DEBUG TEST at Node(k,j,i)=%d,%d,%d: Read { Valid Cell }Ucat[%d][%d][%d]=(%.2f,%.2f,%.2f) and {Out-of-Bounds Cell}Ucat[%d][%d][%d]=(%.2f,%.2f,%.2f)\n",
244 rank, k, j, i,
245 k, j, i_cell_A, ucat_A.x, ucat_A.y, ucat_A.z,
246 k, j, i_cell_B, ucat_B.x, ucat_B.y, ucat_B.z);
247
248 }
249 */
250
251 // Skip processing the unused last node in each dimension.
252 if(i >= IM || j >= JM || k >= KM){
253 continue;
254 }
255 // Loop over the 8 potential cells surrounding node N(k,j,i) and accumulate values.
256 // The index offsets correspond to the relative position of the indices(shifted) that represent cell-centered field values of cells that share the node N(k,j,i) as a corner
257 for (PetscInt dk_offset = -1; dk_offset <= 0; dk_offset++) {
258 for (PetscInt dj_offset = -1; dj_offset <= 0; dj_offset++) {
259 for (PetscInt di_offset = -1; di_offset <= 0; di_offset++) {
260
261 // These are still GLOBAL cell indices
262 PetscInt global_cell_k = k + dk_offset;
263 PetscInt global_cell_j = j + dj_offset;
264 PetscInt global_cell_i = i + di_offset;
265
266 // Check if this corresponds to a valid GLOBAL cell index
267 if (global_cell_i >= 0 && global_cell_i < nCellsX &&
268 global_cell_j >= 0 && global_cell_j < nCellsY &&
269 global_cell_k >= 0 && global_cell_k < nCellsZ)
270 {
271 Cmpnts cell_val = centfield_arr[global_cell_k + 1][global_cell_j + 1][global_cell_i + 1];
272
273 LOG_LOOP_ALLOW_EXACT(LOCAL, LOG_VERBOSE,k,49,"[Rank %d] successful read from [%d][%d][%d] -> (%.2f, %.2f, %.2f)\n",
274 rank,global_cell_k,global_cell_j,global_cell_i,cell_val.x, cell_val.y, cell_val.z);
275
276 sum.x += cell_val.x;
277 sum.y += cell_val.y;
278 sum.z += cell_val.z;
279 count++;
280 }
281 }
282 }
283 }
284
285 PetscInt i_global_write = i; // Global index in GLOBAL array.
286 PetscInt j_global_write = j;
287 PetscInt k_global_write = k;
288
289 // We write directly into the array using the global loop indices.
290 if (count > 0) {
291 corner_arr[k_global_write][j_global_write][i_global_write].x = sum.x / (PetscReal)count;
292 corner_arr[k_global_write][j_global_write][i_global_write].y = sum.y / (PetscReal)count;
293 corner_arr[k_global_write][j_global_write][i_global_write].z = sum.z / (PetscReal)count;
294 } else {
295 // This case should ideally not happen for a valid owned node, but as a failsafe:
296 corner_arr[k_global_write][j_global_write][i_global_write] = (Cmpnts){0.0, 0.0, 0.0};
297 }
298
299 // DEBUG 2
300 /*
301 if(rank == 1){
302 if(i == 11 && j == 11 && k == 49){
303 Cmpnts ucat_node = corner_arr[k][j][i];
304 PetscPrintf(PETSC_COMM_WORLD,"[Rank %d] DEBUG TEST at Node(k,j,i)=%d,%d,%d: Wrote CornerUcat[%d][%d][%d]=(%.2f,%.2f,%.2f)\n",
305 rank, k, j, i,
306 k, j, i, ucat_node.x, ucat_node.y, ucat_node.z);
307 }
308 }
309
310 if(rank == 0 && i == 24 && j == 12 && k == 0){
311 Cmpnts ucat_node = corner_arr[k][j][i];
312 PetscPrintf(PETSC_COMM_WORLD,"[Rank %d] DEBUG TEST at Node(k,j,i)=%d,%d,%d: Wrote CornerUcat[%d][%d][%d]=(%.2f,%.2f,%.2f)\n",
313 rank, k, j, i,
314 k, j, i, ucat_node.x, ucat_node.y, ucat_node.z);
315 }
316 */
317 // LOG_LOOP_ALLOW_EXACT(LOCAL, LOG_VERBOSE,k,48,"[Rank %d] Node(k,j,i)=%d,%d,%d finished loops and write.\n", rank, k, j, i);
318 }
319 }
320 }
322 return 0;
323}
#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 LOG_LOOP_ALLOW_EXACT(scope, level, var, val, fmt,...)
Logs a custom message if a variable equals a specific value.
Definition logging.h:334

◆ InterpolateFieldFromCenterToCorner_Scalar()

PetscErrorCode InterpolateFieldFromCenterToCorner_Scalar ( PetscReal ***  centfield_arr,
PetscReal ***  corner_arr,
UserCtx user 
)

Interpolates a scalar field from cell centers to corner nodes.

This version is adapted to write directly into a ghosted local array obtained from DMDAVecGetArray(), which allows using GLOBAL indices for writing to the OWNED portion of the array.

Parameters
[in]centfield_arrInput: 3D array (ghosted) of scalar data at cell centers, accessed via GLOBAL indices.
[out]corner_arrOutput: 3D array (ghosted) where interpolated node values are stored, also accessed via GLOBAL indices for the owned part.
[in]userUser context containing DMDA information.
Returns
PetscErrorCode 0 on success.

Interpolates a scalar field from cell centers to corner nodes.

Local to this translation unit.

Definition at line 331 of file interpolation.c.

335{
336 PetscErrorCode ierr;
337 DMDALocalInfo info;
338 PetscMPIInt rank;
340 ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
341 ierr = DMDAGetLocalInfo(user->fda, &info); CHKERRQ(ierr);
342
343 // Node ownership range (GLOBAL indices)
344 PetscInt xs_node = info.xs, xm_node = info.xm, xe_node = xs_node + xm_node;
345 PetscInt ys_node = info.ys, ym_node = info.ym, ye_node = ys_node + ym_node;
346 PetscInt zs_node = info.zs, zm_node = info.zm, ze_node = zs_node + zm_node;
347
348 PetscInt nCellsX = info.mx - 2; // Number of cells in x-direction
349 PetscInt nCellsY = info.my - 2; // Number of cells in y-direction
350 PetscInt nCellsZ = info.mz - 2; // Number of cells in z-direction
351
352
353 // Global grid dimensions (used for valid cell check)
354 PetscInt IM = info.mx - 1; // Total nodes in i-direction
355 PetscInt JM = info.my - 1; // Total nodes in j-direction
356 PetscInt KM = info.mz - 1; // Total nodes in k-direction
357
359 "[Rank %d] Starting -- Node ownership k=%d..%d, j=%d..%d, i=%d..%d\n",
360 rank, zs_node, ze_node-1, ys_node, ye_node-1, xs_node, xe_node-1);
361
362 // Loop over the GLOBAL indices of the NODES owned by this processor
363 for (PetscInt k = zs_node; k < ze_node; k++) {
364 for (PetscInt j = ys_node; j < ye_node; j++) {
365 for (PetscInt i = xs_node; i < xe_node; i++) {
366 PetscReal sum = 0.0;
367 PetscInt count = 0;
368
369 // DEBUG 1 TEST
370 /*
371 if(rank == 1 && i == 24 && j == 12 && k == 49){
372 PetscInt i_cell_A = i - 1;
373 PetscInt i_cell_B = i;
374
375 Cmpnts ucat_A = centfield_arr[k][j][i_cell_A]; // 23
376 Cmpnts ucat_B = centfield_arr[k][j][i_cell_B]; // 24 (out-of-bounds if IM=25)
377
378 PetscPrintf(PETSC_COMM_WORLD,"[Rank %d] DEBUG TEST at Node(k,j,i)=%d,%d,%d: Read { Valid Cell }Ucat[%d][%d][%d]=(%.2f,%.2f,%.2f) and {Out-of-Bounds Cell}Ucat[%d][%d][%d]=(%.2f,%.2f,%.2f)\n",
379 rank, k, j, i,
380 k, j, i_cell_A, ucat_A.x, ucat_A.y, ucat_A.z,
381 k, j, i_cell_B, ucat_B.x, ucat_B.y, ucat_B.z);
382
383 }
384 */
385
386 // Skip processing the unused last node in each dimension.
387 if(i >= IM || j >= JM || k >= KM){
388 continue;
389 }
390 // Loop over the 8 potential cells surrounding node N(k,j,i) and accumulate values.
391 // The index offsets correspond to the relative position of the indices(shifted) that represent cell-centered field values of cells that share the node N(k,j,i) as a corner
392 for (PetscInt dk_offset = -1; dk_offset <= 0; dk_offset++) {
393 for (PetscInt dj_offset = -1; dj_offset <= 0; dj_offset++) {
394 for (PetscInt di_offset = -1; di_offset <= 0; di_offset++) {
395
396 // These are still GLOBAL cell indices
397 PetscInt global_cell_k = k + dk_offset;
398 PetscInt global_cell_j = j + dj_offset;
399 PetscInt global_cell_i = i + di_offset;
400
401 // Check if this corresponds to a valid GLOBAL cell index
402 if (global_cell_i >= 0 && global_cell_i < nCellsX &&
403 global_cell_j >= 0 && global_cell_j < nCellsY &&
404 global_cell_k >= 0 && global_cell_k < nCellsZ)
405 {
406 PetscReal cell_val = centfield_arr[global_cell_k + 1][global_cell_j + 1][global_cell_i + 1];
407
408 LOG_LOOP_ALLOW_EXACT(LOCAL, LOG_VERBOSE,k,49,"[Rank %d] successful read from [%d][%d][%d] -> (%.2f)\n",
409 rank,global_cell_k,global_cell_j,global_cell_i,cell_val);
410
411 sum += cell_val;
412 count++;
413 }
414 }
415 }
416 }
417
418 PetscInt i_global_write = i; // Global index in GLOBAL array.
419 PetscInt j_global_write = j;
420 PetscInt k_global_write = k;
421
422 // We write directly into the array using the global loop indices.
423 if (count > 0) {
424 corner_arr[k_global_write][j_global_write][i_global_write] = sum / (PetscReal)count;
425 } else {
426 // This case should ideally not happen for a valid owned node, but as a failsafe:
427 corner_arr[k_global_write][j_global_write][i_global_write] = 0.0;
428 }
429
430 // DEBUG 2
431 /*
432 if(rank == 1){
433 if(i == 11 && j == 11 && k == 49){
434 Cmpnts ucat_node = corner_arr[k][j][i];
435 PetscPrintf(PETSC_COMM_WORLD,"[Rank %d] DEBUG TEST at Node(k,j,i)=%d,%d,%d: Wrote CornerUcat[%d][%d][%d]=(%.2f,%.2f,%.2f)\n",
436 rank, k, j, i,
437 k, j, i, ucat_node.x, ucat_node.y, ucat_node.z);
438 }
439 }
440
441 if(rank == 0 && i == 24 && j == 12 && k == 0){
442 Cmpnts ucat_node = corner_arr[k][j][i];
443 PetscPrintf(PETSC_COMM_WORLD,"[Rank %d] DEBUG TEST at Node(k,j,i)=%d,%d,%d: Wrote CornerUcat[%d][%d][%d]=(%.2f,%.2f,%.2f)\n",
444 rank, k, j, i,
445 k, j, i, ucat_node.x, ucat_node.y, ucat_node.z);
446 }
447 */
448 // LOG_LOOP_ALLOW_EXACT(LOCAL, LOG_VERBOSE,k,48,"[Rank %d] Node(k,j,i)=%d,%d,%d finished loops and write.\n", rank, k, j, i);
449 }
450 }
451 }
453 return 0;
454}

◆ InterpolateCornerToFaceCenter_Scalar()

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.

This routine computes the average of the four corner-node values defining each face of a hexahedral cell:

  • X-faces (perpendicular to X): face between (i-1,i) in X-dir
  • Y-faces (perpendicular to Y): face between (j-1,j) in Y-dir
  • Z-faces (perpendicular to Z): face between (k-1,k) in Z-dir
Parameters
[in]corner_arrGhosted node-centered array (global indexing) from user->fda.
[out]faceX_arrLocal array for X-faces sized [zm][ym][xm+1].
[out]faceY_arrLocal array for Y-faces sized [zm][ym+1][xm].
[out]faceZ_arrLocal array for Z-faces sized [zm+1][ym][xm].
[in]userUser context containing DMDA 'fda' and GetOwnedCellRange.
Returns
PetscErrorCode 0 on success, non-zero on failure.

Interpolates a scalar field from corner nodes to all face centers.

Local to this translation unit.

Definition at line 2069 of file interpolation.c.

2075{
2076 PetscErrorCode ierr;
2077 DMDALocalInfo info;
2078
2079 PetscFunctionBeginUser;
2080
2082
2083 ierr = DMDAGetLocalInfo(user->fda, &info); CHKERRQ(ierr);
2084
2085 // Determine owned-cell ranges based on corner-node ownership
2086 PetscInt xs, xm, ys, ym, zs, zm;
2087 ierr = GetOwnedCellRange(&info, 0, &xs, &xm); CHKERRQ(ierr);
2088 ierr = GetOwnedCellRange(&info, 1, &ys, &ym); CHKERRQ(ierr);
2089 ierr = GetOwnedCellRange(&info, 2, &zs, &zm); CHKERRQ(ierr);
2090
2091 // Global exclusive end indices for cells
2092 PetscInt xe = xs + xm;
2093 PetscInt ye = ys + ym;
2094 PetscInt ze = zs + zm;
2095
2096 // --- X‐faces: loops k=zs..ze-1, j=ys..ye-1, i=xs..xe (xm+1 faces per row) ---
2097 for (PetscInt k = zs; k < ze; ++k) {
2098 PetscInt k_loc = k - zs;
2099 for (PetscInt j = ys; j < ye; ++j) {
2100 PetscInt j_loc = j - ys;
2101 for (PetscInt i = xs; i <= xe; ++i) {
2102 PetscInt i_loc = i - xs; // 0..xm
2103 // Average the four corners of the Y-Z face at X = i
2104 PetscReal sum = corner_arr[k ][j ][i]
2105 + corner_arr[k+1][j ][i]
2106 + corner_arr[k ][j+1][i]
2107 + corner_arr[k+1][j+1][i];
2108 faceX_arr[k_loc][j_loc][i_loc] = sum * 0.25;
2109 }
2110 }
2111 }
2112
2113 // --- Y‐faces: loops k=zs..ze-1, j=ys..ye (ym+1 faces), i=xs..xe-1 ---
2114 for (PetscInt k = zs; k < ze; ++k) {
2115 PetscInt k_loc = k - zs;
2116 for (PetscInt j = ys; j <= ye; ++j) {
2117 PetscInt j_loc = j - ys; // 0..ym
2118 for (PetscInt i = xs; i < xe; ++i) {
2119 PetscInt i_loc = i - xs;
2120 // Average the four corners of the X-Z face at Y = j
2121 PetscReal sum = corner_arr[k ][j][i ]
2122 + corner_arr[k+1][j][i ]
2123 + corner_arr[k ][j][i+1]
2124 + corner_arr[k+1][j][i+1];
2125 faceY_arr[k_loc][j_loc][i_loc] = sum * 0.25;
2126 }
2127 }
2128 }
2129
2130 // --- Z‐faces: loops k=zs..ze (zm+1), j=ys..ye-1, i=xs..xe-1 ---
2131 for (PetscInt k = zs; k <= ze; ++k) {
2132 PetscInt k_loc = k - zs;
2133 for (PetscInt j = ys; j < ye; ++j) {
2134 PetscInt j_loc = j - ys;
2135 for (PetscInt i = xs; i < xe; ++i) {
2136 PetscInt i_loc = i - xs;
2137 // Average the four corners of the X-Y face at Z = k
2138 PetscReal sum = corner_arr[k][j ][i ]
2139 + corner_arr[k][j ][i+1]
2140 + corner_arr[k][j+1][i ]
2141 + corner_arr[k][j+1][i+1];
2142 faceZ_arr[k_loc][j_loc][i_loc] = sum * 0.25;
2143 }
2144 }
2145 }
2146
2148
2149 PetscFunctionReturn(0);
2150}
PetscErrorCode GetOwnedCellRange(const DMDALocalInfo *info_nodes, PetscInt dim, PetscInt *xs_cell_global_out, PetscInt *xm_cell_local_out)
Determines the global starting index and number of CELLS owned by the current processor in a specifie...
Definition setup.c:1883
Here is the call graph for this function:

◆ InterpolateCornerToFaceCenter_Vector()

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.

Identical to the scalar version, except it averages each component of the Cmpnts struct at the four corner-nodes per face.

Parameters
[in]corner_arrGhosted 3-component array (global node indices).
[out]faceX_arrLocal array of Cmpnts for X-faces sized [zm][ym][xm+1].
[out]faceY_arrLocal array of Cmpnts for Y-faces sized [zm][ym+1][xm].
[out]faceZ_arrLocal array of Cmpnts for Z-faces sized [zm+1][ym][xm].
[in]userUser context containing DMDA 'fda'.
Returns
PetscErrorCode 0 on success.

Interpolates a vector field from corner nodes to all face centers.

Local to this translation unit.

Definition at line 2159 of file interpolation.c.

2165{
2166 PetscErrorCode ierr;
2167 DMDALocalInfo info;
2168 PetscMPIInt rank;
2169
2170 PetscFunctionBeginUser;
2171
2173
2174 ierr = MPI_Comm_rank(PETSC_COMM_WORLD, &rank); CHKERRQ(ierr);
2176 "Rank %d starting InterpolateFieldFromCornerToFaceCenter_Vector.\n", rank);
2177
2178 ierr = DMDAGetLocalInfo(user->fda, &info); CHKERRQ(ierr);
2179
2180 PetscInt xs, xm, ys, ym, zs, zm;
2181 ierr = GetOwnedCellRange(&info, 0, &xs, &xm); CHKERRQ(ierr);
2182 ierr = GetOwnedCellRange(&info, 1, &ys, &ym); CHKERRQ(ierr);
2183 ierr = GetOwnedCellRange(&info, 2, &zs, &zm); CHKERRQ(ierr);
2184
2185 PetscInt xe = xs + xm;
2186 PetscInt ye = ys + ym;
2187 PetscInt ze = zs + zm;
2188
2189 // X-faces
2190 for (PetscInt k = zs; k < ze; ++k) {
2191 PetscInt k_loc = k - zs;
2192 for (PetscInt j = ys; j < ye; ++j) {
2193 PetscInt j_loc = j - ys;
2194 for (PetscInt i = xs; i <= xe; ++i) {
2195 PetscInt i_loc = i - xs;
2196 Cmpnts sum = {0,0,0};
2197 sum.x = corner_arr[k ][j ][i].x + corner_arr[k+1][j ][i].x
2198 + corner_arr[k ][j+1][i].x + corner_arr[k+1][j+1][i].x;
2199 sum.y = corner_arr[k ][j ][i].y + corner_arr[k+1][j ][i].y
2200 + corner_arr[k ][j+1][i].y + corner_arr[k+1][j+1][i].y;
2201 sum.z = corner_arr[k ][j ][i].z + corner_arr[k+1][j ][i].z
2202 + corner_arr[k ][j+1][i].z + corner_arr[k+1][j+1][i].z;
2203 faceX_arr[k_loc][j_loc][i_loc].x = sum.x * 0.25;
2204 faceX_arr[k_loc][j_loc][i_loc].y = sum.y * 0.25;
2205 faceX_arr[k_loc][j_loc][i_loc].z = sum.z * 0.25;
2206 }
2207 }
2208 }
2209
2211 "Rank %d x-face Interpolation complete.\n", rank);
2212
2213 // Y-faces
2214 for (PetscInt k = zs; k < ze; ++k) {
2215 PetscInt k_loc = k - zs;
2216 for (PetscInt j = ys; j <= ye; ++j) {
2217 PetscInt j_loc = j - ys;
2218 for (PetscInt i = xs; i < xe; ++i) {
2219 PetscInt i_loc = i - xs;
2220 Cmpnts sum = {0,0,0};
2221 sum.x = corner_arr[k ][j][i ].x + corner_arr[k+1][j][i ].x
2222 + corner_arr[k ][j][i+1].x + corner_arr[k+1][j][i+1].x;
2223 sum.y = corner_arr[k ][j][i ].y + corner_arr[k+1][j][i ].y
2224 + corner_arr[k ][j][i+1].y + corner_arr[k+1][j][i+1].y;
2225 sum.z = corner_arr[k ][j][i ].z + corner_arr[k+1][j][i ].z
2226 + corner_arr[k ][j][i+1].z + corner_arr[k+1][j][i+1].z;
2227 faceY_arr[k_loc][j_loc][i_loc].x = sum.x * 0.25;
2228 faceY_arr[k_loc][j_loc][i_loc].y = sum.y * 0.25;
2229 faceY_arr[k_loc][j_loc][i_loc].z = sum.z * 0.25;
2230 }
2231 }
2232 }
2233
2235 "Rank %d y-face Interpolation complete.\n", rank);
2236
2237 // Z-faces
2238 for (PetscInt k = zs; k <= ze; ++k) {
2239 PetscInt k_loc = k - zs;
2240 for (PetscInt j = ys; j < ye; ++j) {
2241 PetscInt j_loc = j - ys;
2242 for (PetscInt i = xs; i < xe; ++i) {
2243 PetscInt i_loc = i - xs;
2244 Cmpnts sum = {0,0,0};
2245 sum.x = corner_arr[k][j ][i ].x + corner_arr[k][j ][i+1].x
2246 + corner_arr[k][j+1][i ].x + corner_arr[k][j+1][i+1].x;
2247 sum.y = corner_arr[k][j ][i ].y + corner_arr[k][j ][i+1].y
2248 + corner_arr[k][j+1][i ].y + corner_arr[k][j+1][i+1].y;
2249 sum.z = corner_arr[k][j ][i ].z + corner_arr[k][j ][i+1].z
2250 + corner_arr[k][j+1][i ].z + corner_arr[k][j+1][i+1].z;
2251 faceZ_arr[k_loc][j_loc][i_loc].x = sum.x * 0.25;
2252 faceZ_arr[k_loc][j_loc][i_loc].y = sum.y * 0.25;
2253 faceZ_arr[k_loc][j_loc][i_loc].z = sum.z * 0.25;
2254 }
2255 }
2256 }
2257
2259 "Rank %d z-face Interpolation complete.\n", rank);
2260
2262 PetscFunctionReturn(0);
2263}
@ LOG_TRACE
Very fine-grained tracing information for in-depth debugging.
Definition logging.h:32
Here is the call graph for this function: