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

Go to the source code of this file.

Macros

#define __FUNCT__   "AnalyticalSolutionEngine"
 

Functions

PetscErrorCode SetAnalyticalGridInfo (UserCtx *user)
 Sets the grid domain and resolution for analytical solution cases.
 
PetscBool AnalyticalTypeRequiresCustomGeometry (const char *analytical_type)
 Reports whether an analytical type requires custom geometry/decomposition logic.
 
PetscBool AnalyticalTypeSupportsInterpolationError (const char *analytical_type)
 Reports whether an analytical type has a non-trivial velocity field for which interpolation error measurement is meaningful.
 
PetscErrorCode AnalyticalSolutionEngine (SimCtx *simCtx)
 Dispatches to the appropriate analytical solution function based on simulation settings.
 
PetscErrorCode SetAnalyticalSolutionForParticles (Vec tempVec, SimCtx *simCtx)
 Applies the analytical solution to particle velocity vector.
 
PetscErrorCode EvaluateAnalyticalScalarProfile (const SimCtx *simCtx, PetscReal x, PetscReal y, PetscReal z, PetscReal t, PetscReal *value)
 Evaluates the configured verification scalar profile at one physical point.
 
PetscErrorCode SetAnalyticalScalarFieldOnParticles (UserCtx *user, const char *swarm_field_name)
 Writes the configured verification scalar profile onto a particle swarm scalar field.
 
PetscErrorCode SetAnalyticalScalarFieldAtCellCenters (UserCtx *user, Vec targetVec)
 Writes the configured verification scalar profile at physical cell centers into a scalar Vec.
 

Macro Definition Documentation

◆ __FUNCT__

#define __FUNCT__   "AnalyticalSolutionEngine"

Definition at line 73 of file AnalyticalSolutions.h.

Function Documentation

◆ SetAnalyticalGridInfo()

PetscErrorCode SetAnalyticalGridInfo ( UserCtx user)

Sets the grid domain and resolution for analytical solution cases.

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

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

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

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

Sets the grid domain and resolution for analytical solution cases.

Local to this translation unit.

Definition at line 80 of file AnalyticalSolutions.c.

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}
PetscBool AnalyticalTypeRequiresCustomGeometry(const char *analytical_type)
Implementation of AnalyticalTypeRequiresCustomGeometry().
#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_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
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
PetscReal Min_X
Definition variables.h:821
PetscInt KM
Definition variables.h:820
PetscInt _this
Definition variables.h:824
PetscReal ry
Definition variables.h:825
PetscReal Max_Y
Definition variables.h:821
PetscReal rz
Definition variables.h:825
PetscInt JM
Definition variables.h:820
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
PetscInt IM
Definition variables.h:820
PetscReal rx
Definition variables.h:825
PetscReal Max_Z
Definition variables.h:821
The master context for the entire simulation.
Definition variables.h:643
Here is the call graph for this function:
Here is the caller graph for this function:

◆ AnalyticalTypeRequiresCustomGeometry()

PetscBool AnalyticalTypeRequiresCustomGeometry ( const char *  analytical_type)

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

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

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

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

Full API contract (arguments, ownership, side effects) is documented with the header declaration in include/AnalyticalSolutions.h.

See also
AnalyticalTypeRequiresCustomGeometry()

Definition at line 55 of file AnalyticalSolutions.c.

56{
57 if (!analytical_type) return PETSC_FALSE;
58 return (strcmp(analytical_type, "TGV3D") == 0) ? PETSC_TRUE : PETSC_FALSE;
59}
Here is the caller graph for this function:

◆ AnalyticalTypeSupportsInterpolationError()

PetscBool AnalyticalTypeSupportsInterpolationError ( const char *  analytical_type)

Reports whether an analytical type has a non-trivial velocity field for which interpolation error measurement is meaningful.

Types with identically zero velocity (e.g. ZERO_FLOW) return PETSC_FALSE because the interpolation error is trivially zero and uninformative.

Parameters
analytical_typeAnalytical solution type string.
Returns
PETSC_TRUE if interpolation error is meaningful, PETSC_FALSE otherwise.

Reports whether an analytical type has a non-trivial velocity field for which interpolation error measurement is meaningful.

Full API contract (arguments, ownership, side effects) is documented with the header declaration in include/AnalyticalSolutions.h.

See also
AnalyticalTypeSupportsInterpolationError()

Definition at line 67 of file AnalyticalSolutions.c.

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}
Here is the caller graph for this function:

◆ AnalyticalSolutionEngine()

PetscErrorCode AnalyticalSolutionEngine ( SimCtx simCtx)

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

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

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

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

Full API contract (arguments, ownership, side effects) is documented with the header declaration in include/AnalyticalSolutions.h.

See also
AnalyticalSolutionEngine()

Definition at line 184 of file AnalyticalSolutions.c.

185{
186 PetscErrorCode ierr;
187 PetscFunctionBeginUser;
189
190 // -- Before any operation, here is defensive test to ensure that the Corner->Center Interpolation method works
191 //ierr = TestCornerToCenterInterpolation(&(simCtx->usermg.mgctx[simCtx->usermg.mglevels - 1]->user[0]));
192
193 // --- Dispatch based on the string provided by the user ---
194 if (strcmp(simCtx->AnalyticalSolutionType, "TGV3D") == 0) {
195 LOG_ALLOW(GLOBAL, LOG_DEBUG, "Applying Analytical Solution: 3D Taylor-Green Vortex (TGV3D).\n");
196 ierr = SetAnalyticalSolution_TGV3D(simCtx); CHKERRQ(ierr);
197 }
198 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}
static PetscErrorCode SetAnalyticalSolution_ZeroFlow(SimCtx *simCtx)
Internal helper implementation: SetAnalyticalSolution_ZeroFlow().
static PetscErrorCode SetAnalyticalSolution_UniformFlow(SimCtx *simCtx)
Internal helper implementation: SetAnalyticalSolution_UniformFlow().
static PetscErrorCode SetAnalyticalSolution_TGV3D(SimCtx *simCtx)
Internal helper implementation: SetAnalyticalSolution_TGV3D().
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetAnalyticalSolutionForParticles()

PetscErrorCode SetAnalyticalSolutionForParticles ( Vec  tempVec,
SimCtx simCtx 
)

Applies the analytical solution to particle velocity vector.

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

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

Applies the analytical solution to particle velocity vector.

Full API contract (arguments, ownership, side effects) is documented with the header declaration in include/AnalyticalSolutions.h.

See also
SetAnalyticalSolutionForParticles()

Definition at line 564 of file AnalyticalSolutions.c.

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}
static PetscErrorCode SetAnalyticalSolutionForParticles_UniformFlow(Vec tempVec, SimCtx *simCtx)
Internal helper implementation: SetAnalyticalSolutionForParticles_UniformFlow().
static PetscErrorCode SetAnalyticalSolutionForParticles_TGV3D(Vec tempVec, SimCtx *simCtx)
Internal helper implementation: SetAnalyticalSolutionForParticles_TGV3D().
Here is the call graph for this function:
Here is the caller graph for this function:

◆ EvaluateAnalyticalScalarProfile()

PetscErrorCode EvaluateAnalyticalScalarProfile ( const SimCtx simCtx,
PetscReal  x,
PetscReal  y,
PetscReal  z,
PetscReal  t,
PetscReal *  value 
)

Evaluates the configured verification scalar profile at one physical point.

Parameters
[in]simCtxSimulation context providing the scalar verification profile.
[in]xPhysical x coordinate.
[in]yPhysical y coordinate.
[in]zPhysical z coordinate.
[in]tSimulation time.
[out]valueEvaluated scalar value.
Returns
PetscErrorCode 0 on success.

Evaluates the configured verification scalar profile at one physical point.

Full API contract (arguments, ownership, side effects) is documented with the header declaration in include/AnalyticalSolutions.h.

See also
EvaluateAnalyticalScalarProfile()

Definition at line 627 of file AnalyticalSolutions.c.

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}
static PetscErrorCode EvaluateConfiguredScalarProfile(const VerificationScalarConfig *cfg, PetscReal x, PetscReal y, PetscReal z, PetscReal *value)
Internal helper that evaluates the configured scalar verification profile.
VerificationScalarConfig verificationScalar
Definition variables.h:700
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetAnalyticalScalarFieldOnParticles()

PetscErrorCode SetAnalyticalScalarFieldOnParticles ( UserCtx user,
const char *  swarm_field_name 
)

Writes the configured verification scalar profile onto a particle swarm scalar field.

Parameters
[in,out]userBlock-local context providing particle positions.
[in]swarm_field_nameName of the scalar swarm field to overwrite.
Returns
PetscErrorCode 0 on success.

Writes the configured verification scalar profile onto a particle swarm scalar field.

Full API contract (arguments, ownership, side effects) is documented with the header declaration in include/AnalyticalSolutions.h.

See also
SetAnalyticalScalarFieldOnParticles()

Definition at line 653 of file AnalyticalSolutions.c.

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}
PetscErrorCode EvaluateAnalyticalScalarProfile(const SimCtx *simCtx, PetscReal x, PetscReal y, PetscReal z, PetscReal t, PetscReal *value)
Implementation of EvaluateAnalyticalScalarProfile().
PetscReal ti
Definition variables.h:652
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetAnalyticalScalarFieldAtCellCenters()

PetscErrorCode SetAnalyticalScalarFieldAtCellCenters ( UserCtx user,
Vec  targetVec 
)

Writes the configured verification scalar profile at physical cell centers into a scalar Vec.

Parameters
[in]userBlock-local context providing the cell-center coordinates.
[in,out]targetVecScalar Vec on user->da storage to fill with analytical reference values.
Returns
PetscErrorCode 0 on success.

Writes the configured verification scalar profile at physical cell centers into a scalar Vec.

Full API contract (arguments, ownership, side effects) is documented with the header declaration in include/AnalyticalSolutions.h.

See also
SetAnalyticalScalarFieldAtCellCenters()

Definition at line 695 of file AnalyticalSolutions.c.

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}
DMDALocalInfo info
Definition variables.h:818
Vec Cent
Definition variables.h:858
A 3D point or vector with PetscScalar components.
Definition variables.h:100
Here is the call graph for this function:
Here is the caller graph for this function: