Skip to content

adam_prism_cpu_object

ADAM, PRISM (Plasma Research usIng Simulation Methods) equations system class definition, CPU backend.

Source: src/app/prism/cpu/adam_prism_cpu_object.F90

Dependencies

Contents

Variables

NameTypeAttributesDescription
compute_fluxes_maxwellprocedure(compute_convective_fluxes_interface)pointerCompute convective fluxes.

Derived Types

prism_cpu_object

Maxwell equations system class definition, CPU backend.

Inheritance

Extends: prism_common_object

Components

NameTypeAttributesDescription
iotype(io_object)IO handler.
amrtype(amr_object)AMR marker handler.
slicestype(slices_object)Slices handler.
blanesmoantype(blanesmoan_object)Blanes-Moan integrator.
cfmtype(cfm_object)Commutator-Free Magnus integrator.
leapfrogtype(leapfrog_object)Leapfrog integrator.
flailtype(flail_object)Linear algebra methods handler.
wenotype(weno_object)WENO reconstructor.
ibtype(ib_object)Immersed boundary.
rktype(rk_object)Runge-Kutta integrator.
adamtype(adam_object)ADAM, grid + tree + field + maps orchestrator.
fdv_schemecharacter(len=:)allocatableFDV scheme, fd/fv.
fdv_orderinteger(kind=I4P)Order of finite difference/volume schemes, general order.
fdv_half_stencilinteger(kind=I4P)Half stencil length of finite difference/volume schemes.
fdv_half_stencilsinteger(kind=I4P)Half stencil length of fdv schemes for each derivative up to 6.
ngcinteger(kind=I4P)pointerNumber of ghost cells.
niinteger(kind=I4P)pointerNumber of cells in i direction.
njinteger(kind=I4P)pointerNumber of cells in j direction.
nkinteger(kind=I4P)pointerNumber of cells in k direction.
nbinteger(kind=I4P)pointerTotal blocks number for MPI.
blocks_numberinteger(kind=I4P)pointerActual blocks number.
nvinteger(kind=I4P)pointerNumber of variables in q vector.
realm_indexinteger(kind=I4P)Inter-realm seam coupling — currently-published stage buffer key.
stage_activeinteger(kind=I4P)Procedure pointer TBPs for FDV operators (set at initialization by backend).
compute_block_total_variationprocedure(compute_block_total_variation_interface)pass(self), pointerCompute TV.
compute_curlprocedure(compute_curl_interface)pass(self), pointerCompute curl.
compute_derivative1procedure(compute_derivative1_interface)pass(self), pointerCompute deriv1.
compute_derivative2procedure(compute_derivative2_interface)pass(self), pointerCompute deriv2.
compute_derivative4procedure(compute_derivative4_interface)pass(self), pointerCompute deriv4.
compute_divergenceprocedure(compute_divergence_interface)pass(self), pointerCompute dive.
compute_gradientprocedure(compute_gradient_interface)pass(self), pointerCompute grad.
compute_laplacianprocedure(compute_laplacian_interface)pass(self), pointerCompute laplac.
nv_cinteger(kind=I4P)pointerNumber of conservative variables in q vector.
nv_sinteger(kind=I4P)pointerNumber of source variables in q vector.
nv_clinteger(kind=I4P)pointerNumber of divergence cleaning variables in q vector.
qreal(kind=R8P)allocatableConservative cell centered variables.
dqreal(kind=R8P)allocatableResiduals right hand side.
q_picreal(kind=R8P)allocatablePIC variables.
pic_fieldsreal(kind=R8P)allocatableFields value at particle locations.
curlreal(kind=R8P)allocatableCurl fields.
divergencereal(kind=R8P)allocatableDivergence fields.
q_nametype(string)allocatableFields names [1:nv].
dq_nametype(string)allocatableResiduals names [1:nv].
q_pic_nametype(string)allocatablePIC Fields names.
curl_nametype(string)allocatableCurl fields names.
div_nametype(string)allocatableDivergence fields names.
energy_Dreal(kind=R8P)allocatableEnergy of field D, time history.
energy_Breal(kind=R8P)allocatableEnergy of field B, time history.
coil_powerreal(kind=R8P)allocatablePower of coils, time history.
Poynting_fluxreal(kind=R8P)allocatableTotal Poynting flux from boundary, time history.
rms_energy_error_Dreal(kind=R8P)RMS energy error of D field.
rms_energy_error_Breal(kind=R8P)RMS energy error of B field.
max_divergence_Dreal(kind=R8P)Maximum of divergence of D field.
max_divergence_Breal(kind=R8P)Maximum of divergence of B field.
max_divergence_Jreal(kind=R8P)Maximum of divergence of J field.
bctype(prism_bc_object)Boundary conditions.
coiltype(prism_coil_object)Coil source term.
external_fieldstype(prism_external_fields_object)External fields.
fWLayertype(prism_fWLayer_object)Far-field weighted layer.
ictype(prism_ic_object)Initial conditions.
leapfrog_pictype(prism_leapfrog_pic_object)Leapfrog PIC integrator.
numericstype(prism_numerics_object)Numerics configuration.
particle_injectiontype(prism_particle_injection_object)Particle injection.
physicstype(prism_physics_object)Physics configuration.
pictype(prism_pic_object)PIC state.
rk_bctype(prism_rk_bc_object)Runge-Kutta BC handler.
rk_pictype(prism_rk_pic_object)Runge-Kutta PIC handler.
timetype(prism_time_object)Time integration state.
flxyz_creal(kind=R8P)allocatableFluxes at cell center with +/- decomposition for all directions.
flx_freal(kind=R8P)allocatableFluxes along x at cell face.
fly_freal(kind=R8P)allocatableFluxes along y at cell face.
flz_freal(kind=R8P)allocatableFluxes along z at cell face.
compute_residualsprocedure(compute_residuals_interface)pass(self), pointerCompute residuals, space operator.
integrateprocedure(integrate_interface)pass(self), pointerIntegrate, time operator.

Type-Bound Procedures

NameAttributesDescription
load_fdv_from_filepass(self)Load FDV config from file.
finalize_mpi_forestpass(self)Process-global MPI finalize; forest calls it ONCE after all.
after_topology_build_forestpass(self)Backend hook invoked once after the forest builds the seam maps.
close_block_xh5fnopassClose XH5F file block.
close_file_xh5fnopassClose XH5F file.
open_block_xh5fpass(self)Open block file XH5F.
open_file_xh5fpass(self)Open file XH5F.
save_q_xh5fpass(self)Save in XH5F (XDMF/HDF5) format.
allocate_commonpass(self)Allocate common data.
compute_auxiliary_fieldspass(self)Compute auxiliary fields.
initializepass(self)Initialize the equation common data.
load_restart_filespass(self)Load restart files.
save_energy_errorpass(self)Save energy error history.
save_energy_historypass(self)Save energy history.
save_divergence_historypass(self)Save divergence history.
save_restart_filespass(self)Save restart files.
save_xh5fpass(self)Save simulation data in XH5F format.
initialize_coilspass(self)Initialize coils.
set_rectangular_coil_xpass(self)Subroutine to set a rectangular coil source with +-x normal
set_rectangular_coil_ypass(self)Subroutine to set a rectangular coil source with +-y normal
set_rectangular_coil_zpass(self)Subroutine to set a rectangular coil source with +-z normal
set_circular_coil_xpass(self)Subroutine to set a circular coil source with +-x normal
set_circular_coil_ypass(self)Subroutine to set a circular coil source with +-y normal
set_circular_coil_zpass(self)Subroutine to set a circular coil source with +-z normal
set_solenoid_xpass(self)Subroutine to set a solenoid source with +-x normal
set_solenoid_ypass(self)Subroutine to set a solenoid source with +-y normal
set_solenoid_zpass(self)Subroutine to set a solenoid source with +-z normal
set_helicon_coilpass(self)Subroutine to set a helicon coil source.
coupling_descriptor_forestpass(self)Report (scheme_time, rk_scheme, nv) for β admissibility.
amr_updatepass(self)Do AMR update.
mark_by_j_vec_total_variationpass(self)Mark blocks to be refined/derefined by j_vec total variation.
allocate_cpupass(self)Allocate CPU data.
initialize_prismpass(self)Initialize PRSIM equation.
save_residualspass(self)Save residuals history.
save_simulation_datapass(self)Save all simulation data.
apply_fWL_correctionpass(self)Apply fWLayer correction (if present)
compute_coils_currentpass(self)Compute current coils sources.
impose_div_coil_correctionpass(self)Impose coil divergence correction.
set_boundary_conditionspass(self)Set boundary conditions of equation.
compute_residuals_BCpass(self)
update_q_BCpass(self)
set_initial_conditionspass(self)Set initial conditions (and coils) of equation.
update_ghostpass(self)Update ghost cells and set boundary conditions.
compute_field_mean_valuepass(self)Compute field mean value.
initialize_forestpass(self)Invoked by forest%initialize per realm at startup.
compute_local_dt_forestpass(self)Invoked by forest%compute_global_dt during the min reduction.
advance_one_step_forestpass(self)Invoked by forest%evolve_one_step per realm per timestep.
stages_per_step_forestpass(self)Number of integrator stages this realm exposes per step.
open_step_forestpass(self)Per-step prologue (multi-realm path).
begin_stage_forestpass(self)Begin one integrator stage on this realm (multi-realm path).
end_stage_forestpass(self)End the stage: residuals + assignment (multi-realm path).
close_step_forestpass(self)Per-step epilogue (multi-realm path).
post_step_forestpass(self)Invoked by forest%post_step per realm per timestep.
is_done_forestpass(self)Invoked by forest%is_done during the termination reduction.
finalize_forestpass(self)Invoked by forest%finalize per realm at shutdown.
fill_seam_from_peer_forestpass(self)Copy peer's interior into self's ghosts for peer slot p_idx.
apply_reflux_to_stage_forestpass(self)Apply Berger-Colella reflux to self's RK stage buffer.
compute_dtpass(self)Compute time step.
compute_energypass(self)Compute energy.
compute_energy_errorpass(self)Compute energy error.
impose_div_freepass(self)Impose divergence-free property.
impose_ct_correctionpass(self)Impose Constrained Transport correction on q(ivar:ivar+2).
simulatepass(self)Perform the simulation.

Interfaces

compute_residuals_interface

integrate_interface

Subroutines

amr_update

Do AMR update.

fortran
subroutine amr_update(self)

Arguments

NameTypeIntentAttributesDescription
selfclass(prism_cpu_object)inoutThe equation.

Call graph

mark_by_j_vec_total_variation

Mark blocks to be refined/derefined by the value of total variation of j_vec.

fortran
subroutine mark_by_j_vec_total_variation(self, tv_tol, delta_type, delta_fine, delta_coarse, threshold, do_init)

Arguments

NameTypeIntentAttributesDescription
selfclass(prism_cpu_object)inoutThe equation.
tv_tolreal(kind=R8P)inTotal variation tolerance value.
delta_typecharacter(len=*)inDelta criterion type.
delta_finereal(kind=R8P)inMaximum cell delta in fine grids.
delta_coarsereal(kind=R8P)inMinimum cell delta in coarse grids.
thresholdreal(kind=R8P)inoptionalThreshold for sphere proximity.
do_initlogicalinoptionalRe-initialize refinements queries.

Call graph

allocate_cpu

Allocate CPU data.

fortran
subroutine allocate_cpu(self)

Arguments

NameTypeIntentAttributesDescription
selfclass(prism_cpu_object)inoutThe equation.

Call graph

initialize_prism

Initialize PRSIM equation.

fortran
subroutine initialize_prism(self, filename, realms_number)

Arguments

NameTypeIntentAttributesDescription
selfclass(prism_cpu_object)inoutThe equation.
filenamecharacter(len=*)inInput file name.
realms_numberinteger(kind=I4P)inoptionalForest realm count; divides the per-process budget (default 1).

Call graph

save_residuals

Save residuals history.

fortran
subroutine save_residuals(self)

Arguments

NameTypeIntentAttributesDescription
selfclass(prism_cpu_object)inoutThe equation.

Call graph

save_simulation_data

Save all simulation data.

fortran
subroutine save_simulation_data(self)

Arguments

NameTypeIntentAttributesDescription
selfclass(prism_cpu_object)inoutThe equation.

Call graph

compute_coils_current

Compute current coils sources (DC/AC with smooth envelope).

fortran
subroutine compute_coils_current(self, q, gamma)

Arguments

NameTypeIntentAttributesDescription
selfclass(prism_cpu_object)inout
qreal(kind=R8P)inout
gammareal(kind=R8P)inoptional

Call graph

apply_fWL_correction

Apply correction if a fWL is present

fortran
subroutine apply_fWL_correction(self, q)

Arguments

NameTypeIntentAttributesDescription
selfclass(prism_cpu_object)inoutThe equation.
qreal(kind=R8P)inoutConservative variables.

Call graph

set_boundary_conditions

Set boundary conditions of equation.

fortran
subroutine set_boundary_conditions(self, q, s)

Arguments

NameTypeIntentAttributesDescription
selfclass(prism_cpu_object)inoutThe equation.
qreal(kind=R8P)inoutConservative variables.
sinteger(kind=I4P)inoptionalStage counter.

Call graph

compute_residuals_BC

Compute residuals BCs. La sua scrittura si lega all'ordine di interpolazione dell'operatore spaziale. Al momento e' scritto per operatore di secondo ordine (1 punto ghost).

fortran
subroutine compute_residuals_BC(self, s)

Arguments

NameTypeIntentAttributesDescription
selfclass(prism_cpu_object)inoutThe equation.
sinteger(kind=I4P)inStage counter.

Call graph

set_initial_conditions

Set initial conditions and coils on field.

fortran
subroutine set_initial_conditions(self, is_restart)

Arguments

NameTypeIntentAttributesDescription
selfclass(prism_cpu_object)inoutThe equation.
is_restartlogicalinBranching sentinel for restart/non restart path.

Call graph

update_ghost

Update ghost cells (intra-realm only). If not specified all steps are performed, synchronous computation.

Inter-realm seam ghost cells are NO LONGER filled here: the forest's Phase 2 seam fill (between begin_stage and end_stage) populates them via the fill_seam_from_peer_forest TBP operating on the realm's currently-active buffer (selected internally via stage_active). The realm(:) optional dummy that used to thread through this signature has been retired.

fortran
subroutine update_ghost(self, q, step, s)

Arguments

NameTypeIntentAttributesDescription
selfclass(prism_cpu_object)inoutThe equation.
qreal(kind=R8P)inoutConservative variables.
stepinteger(kind=I4P)inoptionalStep to be perfordmed in asyncronous comp.
sinteger(kind=I4P)inoptionalStage counter.

Call graph

initialize_forest

Initialize this realm from scratch: PRISM init, IC injection (or restart load), initial ghost update, initial diagnostics dump, IO files open, plus PIC/leapfrog priming if those schemes are active.

Invoked by forest%initialize. The forest writes self%realm_index = is BEFORE calling this routine, so the body can already read the 1-based forest position from self%realm_index if needed. Pic residual computation

Integration of equations

fortran
subroutine initialize_forest(self, filename, realms_number, memory_avail, nv, verbose)

Arguments

NameTypeIntentAttributesDescription
selfclass(prism_cpu_object)inoutThe realm.
filenamecharacter(len=*)inInput parameters file name.
realms_numberinteger(kind=I4P)inoptionalRealm count; divides the per-process budget.
memory_availreal(kind=R8P)inoptionalPer-process memory budget override.
nvinteger(kind=I4P)inoptionalNumber of field variables override.
verboselogicalinoptionalTrigger verbose output.

Call graph

compute_local_dt_forest

Compute this realm's local stability-limited dt (no MPI reduction).

Invoked by forest%compute_global_dt during the min reduction across all realms in the forest. The reduction itself is the orchestrator's job; this method computes only the value local to self.

fortran
subroutine compute_local_dt_forest(self, dt_local)

Arguments

NameTypeIntentAttributesDescription
selfclass(prism_cpu_object)inThe realm.
dt_localreal(kind=R8P)outLocal stability-limited dt.

Call graph

advance_one_step_forest

Advance this realm by one full timestep of size dt.

Invoked by forest%evolve_one_step once per realm per timestep. Owns the integration itself, i.e. everything that turns q at time t into q at time t + dt.

fortran
subroutine advance_one_step_forest(self, dt)

Arguments

NameTypeIntentAttributesDescription
selfclass(prism_cpu_object)inoutThe realm.
dtreal(kind=R8P)inTimestep size from the forest's global reduction.

Call graph

open_step_forest

Per-step prologue on the multi-realm path: set dt, init RK stages.

Mirrors the head of integrate_rk_ssp: external-field prelude (if active), rk%initialize_stages(q=self%q), and the per-step time bookkeeping advance_one_step_forest does inline (it increment, dt cap for time_max, self%time%dt update). Time advance and progress print run in close_step_forest, mirroring the legacy ordering.

fortran
subroutine open_step_forest(self, dt)

Arguments

NameTypeIntentAttributesDescription
selfclass(prism_cpu_object)inoutThe realm.
dtreal(kind=R8P)inTimestep size from the forest.

Call graph

begin_stage_forest

Begin integrator stage k on the multi-realm path.

Mirrors the FIRST line of integrate_rk_ssp's substage body: rk%compute_stage(s=k, dt=self%time%dt [, phi=ib%phi]). This populates rk%q_rk(:,:,:,:,:,k) from previously computed substages and from self%q, and publishes the stage buffer to peers by setting self%stage_active = k. No ghost reads, no peer-realm access — peer realms may not yet have opened their stage-k buffer when this fires.

fortran
subroutine begin_stage_forest(self, k, K_total, dt, realm)

Arguments

NameTypeIntentAttributesDescription
selfclass(prism_cpu_object)inoutThe realm.
kinteger(kind=I4P)inStage index (1..K_total).
K_totalinteger(kind=I4P)inForest-wide stage count for this step.
dtreal(kind=R8P)inTimestep size from the forest.
realmclass(realm_object)inouttarget, optionalSibling realms (contract parity).

Call graph

end_stage_forest

Close integrator stage k: residuals + stage assignment in one sweep.

fortran
subroutine end_stage_forest(self, k, K_total, dt, realm, flux_register)

Arguments

NameTypeIntentAttributesDescription
selfclass(prism_cpu_object)inoutThe realm.
kinteger(kind=I4P)inStage index (1..K_total).
K_totalinteger(kind=I4P)inForest-wide stage count for this step.
dtreal(kind=R8P)inTimestep size from the forest.
realmclass(realm_object)inouttarget, optionalSibling realms (FNL parity only).
flux_registerclass(flux_register_object)inoutoptionalForest's flux register for FV reflux.

Call graph

close_step_forest

Per-step epilogue on the multi-realm path: q assembly, BC, div-clean, residual save, coil source refresh, time advance, progress print.

Mirrors the post-substage tail of integrate_rk_ssp: rk%update_q + update_q_BC + save_residuals + compute_coils_current + impose_div_free + external_fields add, plus the post-loop time advance + progress print that advance_one_step_forest does inline.

fortran
subroutine close_step_forest(self, dt)

Arguments

NameTypeIntentAttributesDescription
selfclass(prism_cpu_object)inoutThe realm.
dtreal(kind=R8P)inTimestep size from the forest.

Call graph

fill_seam_from_peer_forest

Fill THIS realm's seam ghost cells for peer slot p_idx from peer's interior cells. Single-TBP roundtrip: walks the receiver's sorted seam-map slice for the peer, and for each row copies peer's INTERIOR(b_send, i_send, j_send, k_send, v) → self's GHOST (b_recv, i_recv, j_recv, k_recv, v) on the active conservative buffer (q when stage_active == 0, or rk%q_rk(:,...,stage_active) when active). The PEER's stage_active controls its read source; SELF's stage_active controls its write target.

select type(peer) is confined here — the only place that needs to know the peer's concrete backend type. Cross-backend mixing (CPU realm reading from FNL peer) is not supported and triggers error_stop.

fortran
subroutine fill_seam_from_peer_forest(self, peer, p_idx)

Arguments

NameTypeIntentAttributesDescription
selfclass(prism_cpu_object)inout
peerclass(realm_object)intarget
p_idxinteger(kind=I4P)in

Call graph

apply_reflux_to_stage_forest

PRISM-CPU override of the Berger-Colella reflux correction TBP.

α.r1 cadence: end-of-step barrier. The body fires exactly once per realm per step, at stage == self%rk%nrk (the realm's own final RK substage). Earlier stages return immediately as a no-op. This pairs with M4 (accumulate_seam_fluxes_fv gated on stage_idx == rk%nrk): the realm accumulates its end-of-step face flux into F_coarse(:,:,1) / F_fine_sum(:,:,1) at its final substage, and the same final substage immediately applies the Berger-Colella correction from those accumulators. This is the AMReX Reflux cadence; same-K and asymmetric-K both reduce to the same expression because the accumulators hold the final, committed face flux.

For each face in flux_register where face%coarse_realm == self%realm_index, writes the per-cell correction q_rk(:, seam_cell, b_coarse, stage) += ark(stage) * sign * (dt / dx_coarse) * (F_coarse(:,:,1) - F_fine_sum(:,:,1)) on the matching one-cell-thick seam slice (i ∈ {1, ni} for x-normal faces, etc.). The third-axis index is hardcoded to 1 under α.r1 (register collapsed by M2); the stage parameter is retained on the signature for API stability and used only for the cadence gate and as the write target into self%rk%q_rk(:, ..., stage).

Empty-register / out-of-range guards short-circuit to a no-op, so a realm with no seam faces returns cleanly.

fortran
subroutine apply_reflux_to_stage_forest(self, stage, dt, flux_register)

Arguments

NameTypeIntentAttributesDescription
selfclass(prism_cpu_object)inoutThe realm.
stageinteger(kind=I4P)inIntegrator stage 1..K_total.
dtreal(kind=R8P)inTime step.
flux_registerclass(flux_register_object)inForest's flux register.

Call graph

post_step_forest

Run PRISM-CPU's per-timestep post-step work: state IO, energy diagnostics, divergence diagnostics.

Invoked by forest%post_step. v1 implementation is the verbatim post- step block formerly inline in simulate — every action runs every step, since today's cadence is enforced inside the save_* routines themselves (e.g. save_simulation_data honours io%it_save). The do_* flags are signature-only for now: when the forest takes over cadence the flags will gate the individual calls. For now they are accepted but unused, preserving present-day behavior bit-for-bit.

dt, t, it are not consumed by the current body; they are on the contract so the forest can supply them once it owns time-state (today they are still read from the time module singleton).

realm is accepted on contract parity with FNL; unused on CPU.

fortran
subroutine post_step_forest(self, dt, t, it, do_save_state, do_save_residuals, do_save_restart, do_amr, realm)

Arguments

NameTypeIntentAttributesDescription
selfclass(prism_cpu_object)inoutThe realm.
dtreal(kind=R8P)inTimestep size just advanced.
treal(kind=R8P)inSimulation time after the advance.
itinteger(kind=I4P)inIteration index after the advance.
do_save_statelogicalinoptionalSave state output this step.
do_save_residualslogicalinoptionalSave residuals output this step.
do_save_restartlogicalinoptionalSave restart dump this step.
do_amrlogicalinoptionalRun AMR update this step.
realmclass(realm_object)inouttarget, optionalSibling realms.

Call graph

is_done_forest

Decide whether this realm has reached its local termination criterion.

Invoked by forest%is_done. PRISM-CPU override: matches the legacy condition inline in simulate — terminate when either the simulated time has reached self%time%time_max (time-driven mode, it_max <= 0) or the iteration count has reached self%time%it_max (iteration-driven mode). Today the test reads time-state from the time module singleton, so self is unused; once the forest takes over time bookkeeping the body will consume self%time%* instead.

fortran
subroutine is_done_forest(self, done)

Arguments

NameTypeIntentAttributesDescription
selfclass(prism_cpu_object)inThe realm.
donelogicaloutTrue if this realm is done evolving.

Call graph

finalize_forest

Shut this realm down: final state dump, close residuals/energy/ divergence history files, post-loop divergence diagnostics, finalize MPI handler.

Invoked by forest%finalize. v1 implementation is the verbatim post- loop block formerly inline in simulate. Behavior unchanged.

fortran
subroutine finalize_forest(self)

Arguments

NameTypeIntentAttributesDescription
selfclass(prism_cpu_object)inoutThe realm.

Call graph

compute_dt

Compute the global stability-limited dt and store it on self%time%dt.

Body delegates the local computation to compute_local_dt_forest (orchestrator contract method), then performs the legacy MPI_ALLREDUCE on MPI_COMM_WORLD for backward compatibility with simulate. The forest's compute_global_dt will perform its own reduction, possibly on a per-realm sub-comm; the redundancy disappears once the legacy compute_dt is retired.

fortran
subroutine compute_dt(self)

Arguments

NameTypeIntentAttributesDescription
selfclass(prism_cpu_object)inoutThe equation.

Call graph

compute_energy

Compute energy.

fortran
subroutine compute_energy(self)

Arguments

NameTypeIntentAttributesDescription
selfclass(prism_cpu_object)inoutThe equation.

Call graph

compute_energy_error

Compute energy error.

fortran
subroutine compute_energy_error(self)

Arguments

NameTypeIntentAttributesDescription
selfclass(prism_cpu_object)inoutThe equation.

impose_div_free

Impose divergence-free property.

fortran
subroutine impose_div_free(self)

Arguments

NameTypeIntentAttributesDescription
selfclass(prism_cpu_object)inoutThe equation.

Call graph

impose_ct_correction

Impose Constrained Transport Correction on vectorial variable q(ivar:ivar+2). Note that self%divergence memory is used as buffer, be careful.

fortran
subroutine impose_ct_correction(self, ivar)

Arguments

NameTypeIntentAttributesDescription
selfclass(prism_cpu_object)inoutThe equation.
ivarinteger(kind=I4P)inVariable (start) index in q.

Call graph

simulate

Perform the simulation: legacy single-realm entry point.

fortran
subroutine simulate(self, filename)

Arguments

NameTypeIntentAttributesDescription
selfclass(prism_cpu_object)inoutThe equation.
filenamecharacter(len=*)inInput file name.

Call graph

compute_residuals_fd_centered

Compute residuals of equation, space operator, centered finite difference schemes.

fortran
subroutine compute_residuals_fd_centered(self, q, dq, s, flux_register)

Arguments

NameTypeIntentAttributesDescription
selfclass(prism_cpu_object)inoutThe equation.
qreal(kind=R8P)inoutConservative variables.
dqreal(kind=R8P)inoutResiduals.
sinteger(kind=I4P)inoptionalStage counter.
flux_registerclass(flux_register_object)inoutoptionalFlux register.

Call graph

compute_residuals_fv_centered

Compute residuals of equation, space operator, centered finite volume schemes.

fortran
subroutine compute_residuals_fv_centered(self, q, dq, s, flux_register)

Arguments

NameTypeIntentAttributesDescription
selfclass(prism_cpu_object)inoutThe equation.
qreal(kind=R8P)inoutConservative variables.
dqreal(kind=R8P)inoutResiduals.
sinteger(kind=I4P)inoptionalStage counter.
flux_registerclass(flux_register_object)inoutoptionalForest's flux register; FV reflux.

Call graph

accumulate_seam_fluxes_fv

Accumulate end-of-step face fluxes on inter-realm seam faces of the FV-centered scheme. Gated to the realm's final RK substage under α.r1.

Walks every (block, fec) pair and, for non-zero entries in self%adam%maps%inter_realm_face_register_index(b, fec), packs the appropriate flx_f/fly_f/flz_f face skin into a (nv, nface_cells) slab and routes it to the right register accumulator: positive index → coarse side → accumulate_coarse_flux, negative index → fine side → accumulate_fine_flux. The call site above gates on stage_idx == self%rk%nrk, so this routine fires exactly once per realm per step at the realm's end-of-step.

The register accumulators are sized (nv, nface_cells, 1) under α.r1 (M2 reshape) with the full state-vector width nv (set at register_face time from realm%adam%field%nv). The FV scheme only fills the conservative slice 1:nv_c of the flux arrays; the remaining rows are left at zero in flux_slab so the register's F_coarse - F_fine_sum correction is identically zero outside the conservative variables. The third-axis index passed to accumulate_*_flux is hardcoded to 1 (collapsed register).

For the current same-resolution mirror seam (rmf-2realm), each coarse and fine side carries an identical end-of-step flux at the seam face; F_coarse - F_fine_sum is round-off zero in expectation and the downstream q-correction in forest_object%apply_reflux_corrections is a no-op against the bit-exact regression. The FD-centered case does not reach this routine (different compute_residuals target).

fortran
subroutine accumulate_seam_fluxes_fv(self, ni, nj, nk, nv_c, blocks_number, flux_register)

Arguments

NameTypeIntentAttributesDescription
selfclass(prism_cpu_object)inoutThe realm.
niinteger(kind=I4P)inInterior cell counts.
njinteger(kind=I4P)inInterior cell counts.
nkinteger(kind=I4P)inInterior cell counts.
nv_cinteger(kind=I4P)inNumber of conservative variables.
blocks_numberinteger(kind=I4P)inNumber of local blocks.
flux_registerclass(flux_register_object)inoutForest's flux register.

Call graph

compute_residuals_weno

Compute residuals of equation, space operator, WENO schemes.

fortran
subroutine compute_residuals_weno(self, q, dq, s, flux_register)

Arguments

NameTypeIntentAttributesDescription
selfclass(prism_cpu_object)inoutThe equation.
qreal(kind=R8P)inoutConservative variables.
dqreal(kind=R8P)inoutResiduals.
sinteger(kind=I4P)inoptionalStage counter.
flux_registerclass(flux_register_object)inoutoptionalFlux register.

Call graph

integrate_blanesmoan

Integrate equation, time operator, Yoshida RK scheme.

fortran
subroutine integrate_blanesmoan(self)

Arguments

NameTypeIntentAttributesDescription
selfclass(prism_cpu_object)inoutThe equation.

Call graph

integrate_cfm

Integrate equation, time operator, Commutator-Free Magnus integrator.

fortran
subroutine integrate_cfm(self)

Arguments

NameTypeIntentAttributesDescription
selfclass(prism_cpu_object)inoutThe equation.

Call graph

integrate_leapfrog

Integrate equation, time operator, leapfrog scheme.

fortran
subroutine integrate_leapfrog(self)

Arguments

NameTypeIntentAttributesDescription
selfclass(prism_cpu_object)inoutThe equation.

Call graph

integrate_leapfrog_pic

Integrate equation, time operator, leapfrog scheme for particle in cell Maxwell residuals computation Pic residual computation Integration of equations

fortran
subroutine integrate_leapfrog_pic(self)

Arguments

NameTypeIntentAttributesDescription
selfclass(prism_cpu_object)inoutThe equation.

Call graph

integrate_rk_ls

Integrate equation, time operator, RK classical low storage schemes. Low storage RK working on q_rk(:,:,:,:,:,1)/q as stages, update q in place.

fortran
subroutine integrate_rk_ls(self)

Arguments

NameTypeIntentAttributesDescription
selfclass(prism_cpu_object)inoutThe equation.

Call graph

integrate_rk_ssp

Integrate equation, time operator, SSP RK schemes. SSP RK working on q_rk as stages.

fortran
subroutine integrate_rk_ssp(self)

Arguments

NameTypeIntentAttributesDescription
selfclass(prism_cpu_object)inoutThe equation.

Call graph

integrate_rk_ssp_pic

Integrate equation, time operator, SSP RK schemes. SSP RK working on q_rk as stages.

fortran
subroutine integrate_rk_ssp_pic(self)

Arguments

NameTypeIntentAttributesDescription
selfclass(prism_cpu_object)inoutThe equation.

Call graph

update_q_BC

Update RK q ghost cells.

fortran
subroutine update_q_BC(self, dt, phi)

Arguments

NameTypeIntentAttributesDescription
selfclass(prism_cpu_object)inoutRK object.
dtreal(kind=R8P)inCurrent time step.
phireal(kind=R8P)inoptionalIB distance.

Call graph

integrate_rk_yoshida

Integrate equation, time operator, Yoshida RK scheme.

fortran
subroutine integrate_rk_yoshida(self)

Arguments

NameTypeIntentAttributesDescription
selfclass(prism_cpu_object)inoutThe equation.

Call graph

compute_dxyz_min

Compute minimum dxyz space step.

fortran
subroutine compute_dxyz_min(blocks_number, dxyz, dxyz_min)

Arguments

NameTypeIntentAttributesDescription
blocks_numberinteger(kind=I4P)inNumber of blocks.
dxyzreal(kind=R8P)inXYZ space steps.
dxyz_minreal(kind=R8P)outMinimum space step.

Call graph

compute_fluxes_convective_weno

Compute convective fluxes along direction dir, WENO scheme for space operator.

fortran
subroutine compute_fluxes_convective_weno(dir, blocks_number, ni, nj, nk, ngc, nv_c, weno_s, weno_a, weno_p, weno_d, weno_c, weno_zeps, evmax, erw, elw, chi, q, fluxes)

Arguments

NameTypeIntentAttributesDescription
dirinteger(kind=I4P)inDirection, 1=X, 2=Y, 3=Z.
blocks_numberinteger(kind=I4P)inNumber of blocks.
niinteger(kind=I4P)inGrid cells number in I direction.
njinteger(kind=I4P)inGrid cells number in J direction.
nkinteger(kind=I4P)inGrid cells number in K direction.
ngcinteger(kind=I4P)inGhost cells number.
nv_cinteger(kind=I4P)inNumber of conservative varibales.
weno_sinteger(kind=I4P)inWeno stencils number/dimension.
weno_areal(kind=R8P)inOptimal weights.
weno_preal(kind=R8P)inPolinomials coefficients.
weno_dreal(kind=R8P)inSmoothness indicators coefficients.
weno_creal(kind=R8P)inCentered polinomials coefficients.
weno_zepsreal(kind=R8P)inParameter for avoiding division by zero in computing IS.
evmaxreal(kind=R8P)inMaximum waves speed estimation.
erwreal(kind=R8P)inRight eigenvectors for WENO reconstruction.
elwreal(kind=R8P)inLeft eigenvectors for WENO reconstruction.
chireal(kind=R8P)inSpeed parameter for divergence cleaning.
qreal(kind=R8P)inField variables.
fluxesreal(kind=R8P)inoutFluxes.

Call graph

compute_fluxes_convective_ri_weno

Compute convective fluxes at right interface of b,i,j,k.

fortran
subroutine compute_fluxes_convective_ri_weno(dir, b, i, j, k, ngc, nv_c, weno_s, weno_zeps, weno_a, weno_p, weno_d, weno_c, evmax, erw, elw, chi, si, sir, q, fluxes)

Arguments

NameTypeIntentAttributesDescription
dirinteger(kind=I4P)inDirection, 1=X, 2=Y, 3=Z.
binteger(kind=I4P)inCounter.
iinteger(kind=I4P)inCounter.
jinteger(kind=I4P)inCounter.
kinteger(kind=I4P)inCounter.
ngcinteger(kind=I4P)inGhost cells number.
nv_cinteger(kind=I4P)inNumber of conservative varibales in q vector.
weno_sinteger(kind=I4P)inWeno stencils number/dimension.
weno_zepsreal(kind=R8P)inParameter to avoid division by zero.
weno_areal(kind=R8P)inOptimal weights.
weno_preal(kind=R8P)inPolinomials coefficients.
weno_dreal(kind=R8P)inSmoothness indicators coefficients.
weno_creal(kind=R8P)inCentered polinomials coefficients.
evmaxreal(kind=R8P)inMaximum waves speed estimation.
erwreal(kind=R8P)inRight eigenvectors for WENO reconstruction.
elwreal(kind=R8P)inLeft eigenvectors for WENO reconstruction.
chireal(kind=R8P)inSpeed parameter for divergence cleaning.
siinteger(kind=I4P)inStencil increment.
sirreal(kind=R8P)inStencil increment, real cast.
qreal(kind=R8P)inFields variables.
fluxesreal(kind=R8P)inoutFluxes.

Call graph

compute_fluxes_difference

Compute fluxes difference.

fortran
subroutine compute_fluxes_difference(blocks_number, ni, nj, nk, ngc, nv_c, var_Jx, var_Jy, var_Jz, dx, dy, dz, flx, fly, flz, q, dq)

Arguments

NameTypeIntentAttributesDescription
blocks_numberinteger(kind=I4P)inNumber of blocks.
niinteger(kind=I4P)inGrid cells number in I direction.
njinteger(kind=I4P)inGrid cells number in J direction.
nkinteger(kind=I4P)inGrid cells number in K direction.
ngcinteger(kind=I4P)inGhost cells number.
nv_cinteger(kind=I4P)inNumber of conservative varibales in q.
var_Jxinteger(kind=I4P)inCurrent variable indices.
var_Jyinteger(kind=I4P)inCurrent variable indices.
var_Jzinteger(kind=I4P)inCurrent variable indices.
dxreal(kind=R8P)inSpace steps.
dyreal(kind=R8P)inSpace steps.
dzreal(kind=R8P)inSpace steps.
flxreal(kind=R8P)inX direction fluxes.
flyreal(kind=R8P)inY direction fluxes.
flzreal(kind=R8P)inZ direction fluxes.
qreal(kind=R8P)inFields.
dqreal(kind=R8P)inoutFluxes differences.

Call graph

decompose_fluxes_convective

Decompose convective fluxes.

fortran
subroutine decompose_fluxes_convective(dir, si, sir, b, i, j, k, ngc, nv_c, weno_s, evmax, elw, q, chi, fmpc)

Arguments

NameTypeIntentAttributesDescription
dirinteger(kind=I4P)inDirection, 1=X, 2=Y, 3=Z.
siinteger(kind=I4P)inStencil increment.
sirreal(kind=R8P)inStencil increment, real cast.
binteger(kind=I4P)inCounter.
iinteger(kind=I4P)inCounter.
jinteger(kind=I4P)inCounter.
kinteger(kind=I4P)inCounter.
ngcinteger(kind=I4P)inGhost cells number.
nv_cinteger(kind=I4P)inNumber of conservative varibales in q vector.
weno_sinteger(kind=I4P)inWeno stencils number/dimension.
evmaxreal(kind=R8P)inMaximum eigenvalue.
elwreal(kind=R8P)inLeft eigenvectors for WENO reconstruction.
qreal(kind=R8P)inAuxiliary variables.
chireal(kind=R8P)inSpeed coefficient for D & B div-cleaning.
fmpcreal(kind=R8P)inoutFluxes -+ decomposition in characteristics space.

Call graph

impose_div_coil_correction

Impose Constrained Transport Correction on vectorial variable q(ivar:ivar+2). Note that self%divergence memory is used as buffer, be careful.

fortran
subroutine impose_div_coil_correction(self, ivar, q)

Arguments

NameTypeIntentAttributesDescription
selfclass(prism_cpu_object)inoutThe equation.
ivarinteger(kind=I4P)inVariable (start) index in q.
qreal(kind=R8P)inoutField variables.

Call graph

write_current_behavior_tab

fortran
subroutine write_current_behavior_tab(filename, current_density, time)

Arguments

NameTypeIntentAttributesDescription
filenamecharacter(len=*)in
current_densityreal(kind=R8P)in
timereal(kind=R8P)in

Call graph

write_single_particle_output

fortran
subroutine write_single_particle_output(filename, time, q_pic)

Arguments

NameTypeIntentAttributesDescription
filenamecharacter(len=*)in
timereal(kind=R8P)in
q_picreal(kind=R8P)in

Call graph

compute_field_mean_value

Compute mean value of the field in a certain region of the domain.

fortran
subroutine compute_field_mean_value(self, q, n_x, n_y, n_z, n_b, mean_value)

Arguments

NameTypeIntentAttributesDescription
selfclass(prism_cpu_object)inThe object
qreal(kind=R8P)inField variables.
n_xinteger(kind=I4P)inNumber of cells in each direction and number
n_yinteger(kind=I4P)inNumber of cells in each direction and number
n_zinteger(kind=I4P)inNumber of cells in each direction and number
n_binteger(kind=I4P)inNumber of cells in each direction and number
mean_valuereal(kind=R8P)outMean value of the field out of the fWLayer

Functions

stages_per_step_forest

Number of integrator stages this realm exposes per step.

For the multi-realm path the forest drives the stage loop, so it needs to know K up front. Currently only runge-kutta-ssp-* is split into per-stage TBPs (begin_stage_forest / end_stage_forest); for RK realms K = rk%nrk. Other integrators (Yoshida, Leapfrog, Blanes-Moan, CFM) will error-stop on the multi-realm path until they are split as well.

Returns: integer(kind=I4P)

fortran
function stages_per_step_forest(self) result(K)

Arguments

NameTypeIntentAttributesDescription
selfclass(prism_cpu_object)inThe realm.

Call graph

dotproduct

Compute the scalar (dot) product.

Returns: real(kind=R8P)

fortran
function dotproduct(a, b) result(dot)

Arguments

NameTypeIntentAttributesDescription
areal(kind=R8P)inLeft hand side.
breal(kind=R8P)inLeft hand side.

Call graph

crossproduct

Returns: real(kind=R8P)

fortran
function crossproduct(a, b) result(cross)

Arguments

NameTypeIntentAttributesDescription
areal(kind=R8P)inLeft hand side.
breal(kind=R8P)inLeft hand side.

Call graph