Appearance
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
- prism_cpu_object
- compute_residuals_interface
- integrate_interface
- amr_update
- mark_by_j_vec_total_variation
- allocate_cpu
- initialize_prism
- save_residuals
- save_simulation_data
- compute_coils_current
- apply_fWL_correction
- set_boundary_conditions
- compute_residuals_BC
- set_initial_conditions
- update_ghost
- initialize_forest
- compute_local_dt_forest
- advance_one_step_forest
- open_step_forest
- begin_stage_forest
- end_stage_forest
- close_step_forest
- fill_seam_from_peer_forest
- apply_reflux_to_stage_forest
- post_step_forest
- is_done_forest
- finalize_forest
- compute_dt
- compute_energy
- compute_energy_error
- impose_div_free
- impose_ct_correction
- simulate
- compute_residuals_fd_centered
- compute_residuals_fv_centered
- accumulate_seam_fluxes_fv
- compute_residuals_weno
- integrate_blanesmoan
- integrate_cfm
- integrate_leapfrog
- integrate_leapfrog_pic
- integrate_rk_ls
- integrate_rk_ssp
- integrate_rk_ssp_pic
- update_q_BC
- integrate_rk_yoshida
- compute_dxyz_min
- compute_fluxes_convective_weno
- compute_fluxes_convective_ri_weno
- compute_fluxes_difference
- decompose_fluxes_convective
- impose_div_coil_correction
- write_current_behavior_tab
- write_single_particle_output
- compute_field_mean_value
- stages_per_step_forest
- dotproduct
- crossproduct
Variables
| Name | Type | Attributes | Description |
|---|---|---|---|
compute_fluxes_maxwell | procedure(compute_convective_fluxes_interface) | pointer | Compute convective fluxes. |
Derived Types
prism_cpu_object
Maxwell equations system class definition, CPU backend.
Inheritance
Extends: prism_common_object
Components
| Name | Type | Attributes | Description |
|---|---|---|---|
io | type(io_object) | IO handler. | |
amr | type(amr_object) | AMR marker handler. | |
slices | type(slices_object) | Slices handler. | |
blanesmoan | type(blanesmoan_object) | Blanes-Moan integrator. | |
cfm | type(cfm_object) | Commutator-Free Magnus integrator. | |
leapfrog | type(leapfrog_object) | Leapfrog integrator. | |
flail | type(flail_object) | Linear algebra methods handler. | |
weno | type(weno_object) | WENO reconstructor. | |
ib | type(ib_object) | Immersed boundary. | |
rk | type(rk_object) | Runge-Kutta integrator. | |
adam | type(adam_object) | ADAM, grid + tree + field + maps orchestrator. | |
fdv_scheme | character(len=:) | allocatable | FDV scheme, fd/fv. |
fdv_order | integer(kind=I4P) | Order of finite difference/volume schemes, general order. | |
fdv_half_stencil | integer(kind=I4P) | Half stencil length of finite difference/volume schemes. | |
fdv_half_stencils | integer(kind=I4P) | Half stencil length of fdv schemes for each derivative up to 6. | |
ngc | integer(kind=I4P) | pointer | Number of ghost cells. |
ni | integer(kind=I4P) | pointer | Number of cells in i direction. |
nj | integer(kind=I4P) | pointer | Number of cells in j direction. |
nk | integer(kind=I4P) | pointer | Number of cells in k direction. |
nb | integer(kind=I4P) | pointer | Total blocks number for MPI. |
blocks_number | integer(kind=I4P) | pointer | Actual blocks number. |
nv | integer(kind=I4P) | pointer | Number of variables in q vector. |
realm_index | integer(kind=I4P) | Inter-realm seam coupling — currently-published stage buffer key. | |
stage_active | integer(kind=I4P) | Procedure pointer TBPs for FDV operators (set at initialization by backend). | |
compute_block_total_variation | procedure(compute_block_total_variation_interface) | pass(self), pointer | Compute TV. |
compute_curl | procedure(compute_curl_interface) | pass(self), pointer | Compute curl. |
compute_derivative1 | procedure(compute_derivative1_interface) | pass(self), pointer | Compute deriv1. |
compute_derivative2 | procedure(compute_derivative2_interface) | pass(self), pointer | Compute deriv2. |
compute_derivative4 | procedure(compute_derivative4_interface) | pass(self), pointer | Compute deriv4. |
compute_divergence | procedure(compute_divergence_interface) | pass(self), pointer | Compute dive. |
compute_gradient | procedure(compute_gradient_interface) | pass(self), pointer | Compute grad. |
compute_laplacian | procedure(compute_laplacian_interface) | pass(self), pointer | Compute laplac. |
nv_c | integer(kind=I4P) | pointer | Number of conservative variables in q vector. |
nv_s | integer(kind=I4P) | pointer | Number of source variables in q vector. |
nv_cl | integer(kind=I4P) | pointer | Number of divergence cleaning variables in q vector. |
q | real(kind=R8P) | allocatable | Conservative cell centered variables. |
dq | real(kind=R8P) | allocatable | Residuals right hand side. |
q_pic | real(kind=R8P) | allocatable | PIC variables. |
pic_fields | real(kind=R8P) | allocatable | Fields value at particle locations. |
curl | real(kind=R8P) | allocatable | Curl fields. |
divergence | real(kind=R8P) | allocatable | Divergence fields. |
q_name | type(string) | allocatable | Fields names [1:nv]. |
dq_name | type(string) | allocatable | Residuals names [1:nv]. |
q_pic_name | type(string) | allocatable | PIC Fields names. |
curl_name | type(string) | allocatable | Curl fields names. |
div_name | type(string) | allocatable | Divergence fields names. |
energy_D | real(kind=R8P) | allocatable | Energy of field D, time history. |
energy_B | real(kind=R8P) | allocatable | Energy of field B, time history. |
coil_power | real(kind=R8P) | allocatable | Power of coils, time history. |
Poynting_flux | real(kind=R8P) | allocatable | Total Poynting flux from boundary, time history. |
rms_energy_error_D | real(kind=R8P) | RMS energy error of D field. | |
rms_energy_error_B | real(kind=R8P) | RMS energy error of B field. | |
max_divergence_D | real(kind=R8P) | Maximum of divergence of D field. | |
max_divergence_B | real(kind=R8P) | Maximum of divergence of B field. | |
max_divergence_J | real(kind=R8P) | Maximum of divergence of J field. | |
bc | type(prism_bc_object) | Boundary conditions. | |
coil | type(prism_coil_object) | Coil source term. | |
external_fields | type(prism_external_fields_object) | External fields. | |
fWLayer | type(prism_fWLayer_object) | Far-field weighted layer. | |
ic | type(prism_ic_object) | Initial conditions. | |
leapfrog_pic | type(prism_leapfrog_pic_object) | Leapfrog PIC integrator. | |
numerics | type(prism_numerics_object) | Numerics configuration. | |
particle_injection | type(prism_particle_injection_object) | Particle injection. | |
physics | type(prism_physics_object) | Physics configuration. | |
pic | type(prism_pic_object) | PIC state. | |
rk_bc | type(prism_rk_bc_object) | Runge-Kutta BC handler. | |
rk_pic | type(prism_rk_pic_object) | Runge-Kutta PIC handler. | |
time | type(prism_time_object) | Time integration state. | |
flxyz_c | real(kind=R8P) | allocatable | Fluxes at cell center with +/- decomposition for all directions. |
flx_f | real(kind=R8P) | allocatable | Fluxes along x at cell face. |
fly_f | real(kind=R8P) | allocatable | Fluxes along y at cell face. |
flz_f | real(kind=R8P) | allocatable | Fluxes along z at cell face. |
compute_residuals | procedure(compute_residuals_interface) | pass(self), pointer | Compute residuals, space operator. |
integrate | procedure(integrate_interface) | pass(self), pointer | Integrate, time operator. |
Type-Bound Procedures
| Name | Attributes | Description |
|---|---|---|
load_fdv_from_file | pass(self) | Load FDV config from file. |
finalize_mpi_forest | pass(self) | Process-global MPI finalize; forest calls it ONCE after all. |
after_topology_build_forest | pass(self) | Backend hook invoked once after the forest builds the seam maps. |
close_block_xh5f | nopass | Close XH5F file block. |
close_file_xh5f | nopass | Close XH5F file. |
open_block_xh5f | pass(self) | Open block file XH5F. |
open_file_xh5f | pass(self) | Open file XH5F. |
save_q_xh5f | pass(self) | Save in XH5F (XDMF/HDF5) format. |
allocate_common | pass(self) | Allocate common data. |
compute_auxiliary_fields | pass(self) | Compute auxiliary fields. |
initialize | pass(self) | Initialize the equation common data. |
load_restart_files | pass(self) | Load restart files. |
save_energy_error | pass(self) | Save energy error history. |
save_energy_history | pass(self) | Save energy history. |
save_divergence_history | pass(self) | Save divergence history. |
save_restart_files | pass(self) | Save restart files. |
save_xh5f | pass(self) | Save simulation data in XH5F format. |
initialize_coils | pass(self) | Initialize coils. |
set_rectangular_coil_x | pass(self) | Subroutine to set a rectangular coil source with +-x normal |
set_rectangular_coil_y | pass(self) | Subroutine to set a rectangular coil source with +-y normal |
set_rectangular_coil_z | pass(self) | Subroutine to set a rectangular coil source with +-z normal |
set_circular_coil_x | pass(self) | Subroutine to set a circular coil source with +-x normal |
set_circular_coil_y | pass(self) | Subroutine to set a circular coil source with +-y normal |
set_circular_coil_z | pass(self) | Subroutine to set a circular coil source with +-z normal |
set_solenoid_x | pass(self) | Subroutine to set a solenoid source with +-x normal |
set_solenoid_y | pass(self) | Subroutine to set a solenoid source with +-y normal |
set_solenoid_z | pass(self) | Subroutine to set a solenoid source with +-z normal |
set_helicon_coil | pass(self) | Subroutine to set a helicon coil source. |
coupling_descriptor_forest | pass(self) | Report (scheme_time, rk_scheme, nv) for β admissibility. |
amr_update | pass(self) | Do AMR update. |
mark_by_j_vec_total_variation | pass(self) | Mark blocks to be refined/derefined by j_vec total variation. |
allocate_cpu | pass(self) | Allocate CPU data. |
initialize_prism | pass(self) | Initialize PRSIM equation. |
save_residuals | pass(self) | Save residuals history. |
save_simulation_data | pass(self) | Save all simulation data. |
apply_fWL_correction | pass(self) | Apply fWLayer correction (if present) |
compute_coils_current | pass(self) | Compute current coils sources. |
impose_div_coil_correction | pass(self) | Impose coil divergence correction. |
set_boundary_conditions | pass(self) | Set boundary conditions of equation. |
compute_residuals_BC | pass(self) | |
update_q_BC | pass(self) | |
set_initial_conditions | pass(self) | Set initial conditions (and coils) of equation. |
update_ghost | pass(self) | Update ghost cells and set boundary conditions. |
compute_field_mean_value | pass(self) | Compute field mean value. |
initialize_forest | pass(self) | Invoked by forest%initialize per realm at startup. |
compute_local_dt_forest | pass(self) | Invoked by forest%compute_global_dt during the min reduction. |
advance_one_step_forest | pass(self) | Invoked by forest%evolve_one_step per realm per timestep. |
stages_per_step_forest | pass(self) | Number of integrator stages this realm exposes per step. |
open_step_forest | pass(self) | Per-step prologue (multi-realm path). |
begin_stage_forest | pass(self) | Begin one integrator stage on this realm (multi-realm path). |
end_stage_forest | pass(self) | End the stage: residuals + assignment (multi-realm path). |
close_step_forest | pass(self) | Per-step epilogue (multi-realm path). |
post_step_forest | pass(self) | Invoked by forest%post_step per realm per timestep. |
is_done_forest | pass(self) | Invoked by forest%is_done during the termination reduction. |
finalize_forest | pass(self) | Invoked by forest%finalize per realm at shutdown. |
fill_seam_from_peer_forest | pass(self) | Copy peer's interior into self's ghosts for peer slot p_idx. |
apply_reflux_to_stage_forest | pass(self) | Apply Berger-Colella reflux to self's RK stage buffer. |
compute_dt | pass(self) | Compute time step. |
compute_energy | pass(self) | Compute energy. |
compute_energy_error | pass(self) | Compute energy error. |
impose_div_free | pass(self) | Impose divergence-free property. |
impose_ct_correction | pass(self) | Impose Constrained Transport correction on q(ivar:ivar+2). |
simulate | pass(self) | Perform the simulation. |
Interfaces
compute_residuals_interface
integrate_interface
Subroutines
amr_update
Do AMR update.
fortran
subroutine amr_update(self)Arguments
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
self | class(prism_cpu_object) | inout | The 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
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
self | class(prism_cpu_object) | inout | The equation. | |
tv_tol | real(kind=R8P) | in | Total variation tolerance value. | |
delta_type | character(len=*) | in | Delta criterion type. | |
delta_fine | real(kind=R8P) | in | Maximum cell delta in fine grids. | |
delta_coarse | real(kind=R8P) | in | Minimum cell delta in coarse grids. | |
threshold | real(kind=R8P) | in | optional | Threshold for sphere proximity. |
do_init | logical | in | optional | Re-initialize refinements queries. |
Call graph
allocate_cpu
Allocate CPU data.
fortran
subroutine allocate_cpu(self)Arguments
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
self | class(prism_cpu_object) | inout | The equation. |
Call graph
initialize_prism
Initialize PRSIM equation.
fortran
subroutine initialize_prism(self, filename, realms_number)Arguments
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
self | class(prism_cpu_object) | inout | The equation. | |
filename | character(len=*) | in | Input file name. | |
realms_number | integer(kind=I4P) | in | optional | Forest realm count; divides the per-process budget (default 1). |
Call graph
save_residuals
Save residuals history.
fortran
subroutine save_residuals(self)Arguments
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
self | class(prism_cpu_object) | inout | The equation. |
Call graph
save_simulation_data
Save all simulation data.
fortran
subroutine save_simulation_data(self)Arguments
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
self | class(prism_cpu_object) | inout | The equation. |
Call graph
compute_coils_current
Compute current coils sources (DC/AC with smooth envelope).
fortran
subroutine compute_coils_current(self, q, gamma)Arguments
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
self | class(prism_cpu_object) | inout | ||
q | real(kind=R8P) | inout | ||
gamma | real(kind=R8P) | in | optional |
Call graph
apply_fWL_correction
Apply correction if a fWL is present
fortran
subroutine apply_fWL_correction(self, q)Arguments
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
self | class(prism_cpu_object) | inout | The equation. | |
q | real(kind=R8P) | inout | Conservative variables. |
Call graph
set_boundary_conditions
Set boundary conditions of equation.
fortran
subroutine set_boundary_conditions(self, q, s)Arguments
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
self | class(prism_cpu_object) | inout | The equation. | |
q | real(kind=R8P) | inout | Conservative variables. | |
s | integer(kind=I4P) | in | optional | Stage 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
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
self | class(prism_cpu_object) | inout | The equation. | |
s | integer(kind=I4P) | in | Stage counter. |
Call graph
set_initial_conditions
Set initial conditions and coils on field.
fortran
subroutine set_initial_conditions(self, is_restart)Arguments
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
self | class(prism_cpu_object) | inout | The equation. | |
is_restart | logical | in | Branching 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
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
self | class(prism_cpu_object) | inout | The equation. | |
q | real(kind=R8P) | inout | Conservative variables. | |
step | integer(kind=I4P) | in | optional | Step to be perfordmed in asyncronous comp. |
s | integer(kind=I4P) | in | optional | Stage 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
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
self | class(prism_cpu_object) | inout | The realm. | |
filename | character(len=*) | in | Input parameters file name. | |
realms_number | integer(kind=I4P) | in | optional | Realm count; divides the per-process budget. |
memory_avail | real(kind=R8P) | in | optional | Per-process memory budget override. |
nv | integer(kind=I4P) | in | optional | Number of field variables override. |
verbose | logical | in | optional | Trigger 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
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
self | class(prism_cpu_object) | in | The realm. | |
dt_local | real(kind=R8P) | out | Local 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
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
self | class(prism_cpu_object) | inout | The realm. | |
dt | real(kind=R8P) | in | Timestep 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
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
self | class(prism_cpu_object) | inout | The realm. | |
dt | real(kind=R8P) | in | Timestep 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
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
self | class(prism_cpu_object) | inout | The realm. | |
k | integer(kind=I4P) | in | Stage index (1..K_total). | |
K_total | integer(kind=I4P) | in | Forest-wide stage count for this step. | |
dt | real(kind=R8P) | in | Timestep size from the forest. | |
realm | class(realm_object) | inout | target, optional | Sibling 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
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
self | class(prism_cpu_object) | inout | The realm. | |
k | integer(kind=I4P) | in | Stage index (1..K_total). | |
K_total | integer(kind=I4P) | in | Forest-wide stage count for this step. | |
dt | real(kind=R8P) | in | Timestep size from the forest. | |
realm | class(realm_object) | inout | target, optional | Sibling realms (FNL parity only). |
flux_register | class(flux_register_object) | inout | optional | Forest'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
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
self | class(prism_cpu_object) | inout | The realm. | |
dt | real(kind=R8P) | in | Timestep 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
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
self | class(prism_cpu_object) | inout | ||
peer | class(realm_object) | in | target | |
p_idx | integer(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
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
self | class(prism_cpu_object) | inout | The realm. | |
stage | integer(kind=I4P) | in | Integrator stage 1..K_total. | |
dt | real(kind=R8P) | in | Time step. | |
flux_register | class(flux_register_object) | in | Forest'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
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
self | class(prism_cpu_object) | inout | The realm. | |
dt | real(kind=R8P) | in | Timestep size just advanced. | |
t | real(kind=R8P) | in | Simulation time after the advance. | |
it | integer(kind=I4P) | in | Iteration index after the advance. | |
do_save_state | logical | in | optional | Save state output this step. |
do_save_residuals | logical | in | optional | Save residuals output this step. |
do_save_restart | logical | in | optional | Save restart dump this step. |
do_amr | logical | in | optional | Run AMR update this step. |
realm | class(realm_object) | inout | target, optional | Sibling 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
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
self | class(prism_cpu_object) | in | The realm. | |
done | logical | out | True 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
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
self | class(prism_cpu_object) | inout | The 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
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
self | class(prism_cpu_object) | inout | The equation. |
Call graph
compute_energy
Compute energy.
fortran
subroutine compute_energy(self)Arguments
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
self | class(prism_cpu_object) | inout | The equation. |
Call graph
compute_energy_error
Compute energy error.
fortran
subroutine compute_energy_error(self)Arguments
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
self | class(prism_cpu_object) | inout | The equation. |
impose_div_free
Impose divergence-free property.
fortran
subroutine impose_div_free(self)Arguments
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
self | class(prism_cpu_object) | inout | The 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
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
self | class(prism_cpu_object) | inout | The equation. | |
ivar | integer(kind=I4P) | in | Variable (start) index in q. |
Call graph
simulate
Perform the simulation: legacy single-realm entry point.
fortran
subroutine simulate(self, filename)Arguments
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
self | class(prism_cpu_object) | inout | The equation. | |
filename | character(len=*) | in | Input 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
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
self | class(prism_cpu_object) | inout | The equation. | |
q | real(kind=R8P) | inout | Conservative variables. | |
dq | real(kind=R8P) | inout | Residuals. | |
s | integer(kind=I4P) | in | optional | Stage counter. |
flux_register | class(flux_register_object) | inout | optional | Flux 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
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
self | class(prism_cpu_object) | inout | The equation. | |
q | real(kind=R8P) | inout | Conservative variables. | |
dq | real(kind=R8P) | inout | Residuals. | |
s | integer(kind=I4P) | in | optional | Stage counter. |
flux_register | class(flux_register_object) | inout | optional | Forest'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
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
self | class(prism_cpu_object) | inout | The realm. | |
ni | integer(kind=I4P) | in | Interior cell counts. | |
nj | integer(kind=I4P) | in | Interior cell counts. | |
nk | integer(kind=I4P) | in | Interior cell counts. | |
nv_c | integer(kind=I4P) | in | Number of conservative variables. | |
blocks_number | integer(kind=I4P) | in | Number of local blocks. | |
flux_register | class(flux_register_object) | inout | Forest'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
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
self | class(prism_cpu_object) | inout | The equation. | |
q | real(kind=R8P) | inout | Conservative variables. | |
dq | real(kind=R8P) | inout | Residuals. | |
s | integer(kind=I4P) | in | optional | Stage counter. |
flux_register | class(flux_register_object) | inout | optional | Flux register. |
Call graph
integrate_blanesmoan
Integrate equation, time operator, Yoshida RK scheme.
fortran
subroutine integrate_blanesmoan(self)Arguments
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
self | class(prism_cpu_object) | inout | The equation. |
Call graph
integrate_cfm
Integrate equation, time operator, Commutator-Free Magnus integrator.
fortran
subroutine integrate_cfm(self)Arguments
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
self | class(prism_cpu_object) | inout | The equation. |
Call graph
integrate_leapfrog
Integrate equation, time operator, leapfrog scheme.
fortran
subroutine integrate_leapfrog(self)Arguments
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
self | class(prism_cpu_object) | inout | The 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
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
self | class(prism_cpu_object) | inout | The 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
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
self | class(prism_cpu_object) | inout | The 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
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
self | class(prism_cpu_object) | inout | The 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
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
self | class(prism_cpu_object) | inout | The equation. |
Call graph
update_q_BC
Update RK q ghost cells.
fortran
subroutine update_q_BC(self, dt, phi)Arguments
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
self | class(prism_cpu_object) | inout | RK object. | |
dt | real(kind=R8P) | in | Current time step. | |
phi | real(kind=R8P) | in | optional | IB distance. |
Call graph
integrate_rk_yoshida
Integrate equation, time operator, Yoshida RK scheme.
fortran
subroutine integrate_rk_yoshida(self)Arguments
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
self | class(prism_cpu_object) | inout | The equation. |
Call graph
compute_dxyz_min
Compute minimum dxyz space step.
fortran
subroutine compute_dxyz_min(blocks_number, dxyz, dxyz_min)Arguments
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
blocks_number | integer(kind=I4P) | in | Number of blocks. | |
dxyz | real(kind=R8P) | in | XYZ space steps. | |
dxyz_min | real(kind=R8P) | out | Minimum 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
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
dir | integer(kind=I4P) | in | Direction, 1=X, 2=Y, 3=Z. | |
blocks_number | integer(kind=I4P) | in | Number of blocks. | |
ni | integer(kind=I4P) | in | Grid cells number in I direction. | |
nj | integer(kind=I4P) | in | Grid cells number in J direction. | |
nk | integer(kind=I4P) | in | Grid cells number in K direction. | |
ngc | integer(kind=I4P) | in | Ghost cells number. | |
nv_c | integer(kind=I4P) | in | Number of conservative varibales. | |
weno_s | integer(kind=I4P) | in | Weno stencils number/dimension. | |
weno_a | real(kind=R8P) | in | Optimal weights. | |
weno_p | real(kind=R8P) | in | Polinomials coefficients. | |
weno_d | real(kind=R8P) | in | Smoothness indicators coefficients. | |
weno_c | real(kind=R8P) | in | Centered polinomials coefficients. | |
weno_zeps | real(kind=R8P) | in | Parameter for avoiding division by zero in computing IS. | |
evmax | real(kind=R8P) | in | Maximum waves speed estimation. | |
erw | real(kind=R8P) | in | Right eigenvectors for WENO reconstruction. | |
elw | real(kind=R8P) | in | Left eigenvectors for WENO reconstruction. | |
chi | real(kind=R8P) | in | Speed parameter for divergence cleaning. | |
q | real(kind=R8P) | in | Field variables. | |
fluxes | real(kind=R8P) | inout | Fluxes. |
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
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
dir | integer(kind=I4P) | in | Direction, 1=X, 2=Y, 3=Z. | |
b | integer(kind=I4P) | in | Counter. | |
i | integer(kind=I4P) | in | Counter. | |
j | integer(kind=I4P) | in | Counter. | |
k | integer(kind=I4P) | in | Counter. | |
ngc | integer(kind=I4P) | in | Ghost cells number. | |
nv_c | integer(kind=I4P) | in | Number of conservative varibales in q vector. | |
weno_s | integer(kind=I4P) | in | Weno stencils number/dimension. | |
weno_zeps | real(kind=R8P) | in | Parameter to avoid division by zero. | |
weno_a | real(kind=R8P) | in | Optimal weights. | |
weno_p | real(kind=R8P) | in | Polinomials coefficients. | |
weno_d | real(kind=R8P) | in | Smoothness indicators coefficients. | |
weno_c | real(kind=R8P) | in | Centered polinomials coefficients. | |
evmax | real(kind=R8P) | in | Maximum waves speed estimation. | |
erw | real(kind=R8P) | in | Right eigenvectors for WENO reconstruction. | |
elw | real(kind=R8P) | in | Left eigenvectors for WENO reconstruction. | |
chi | real(kind=R8P) | in | Speed parameter for divergence cleaning. | |
si | integer(kind=I4P) | in | Stencil increment. | |
sir | real(kind=R8P) | in | Stencil increment, real cast. | |
q | real(kind=R8P) | in | Fields variables. | |
fluxes | real(kind=R8P) | inout | Fluxes. |
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
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
blocks_number | integer(kind=I4P) | in | Number of blocks. | |
ni | integer(kind=I4P) | in | Grid cells number in I direction. | |
nj | integer(kind=I4P) | in | Grid cells number in J direction. | |
nk | integer(kind=I4P) | in | Grid cells number in K direction. | |
ngc | integer(kind=I4P) | in | Ghost cells number. | |
nv_c | integer(kind=I4P) | in | Number of conservative varibales in q. | |
var_Jx | integer(kind=I4P) | in | Current variable indices. | |
var_Jy | integer(kind=I4P) | in | Current variable indices. | |
var_Jz | integer(kind=I4P) | in | Current variable indices. | |
dx | real(kind=R8P) | in | Space steps. | |
dy | real(kind=R8P) | in | Space steps. | |
dz | real(kind=R8P) | in | Space steps. | |
flx | real(kind=R8P) | in | X direction fluxes. | |
fly | real(kind=R8P) | in | Y direction fluxes. | |
flz | real(kind=R8P) | in | Z direction fluxes. | |
q | real(kind=R8P) | in | Fields. | |
dq | real(kind=R8P) | inout | Fluxes 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
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
dir | integer(kind=I4P) | in | Direction, 1=X, 2=Y, 3=Z. | |
si | integer(kind=I4P) | in | Stencil increment. | |
sir | real(kind=R8P) | in | Stencil increment, real cast. | |
b | integer(kind=I4P) | in | Counter. | |
i | integer(kind=I4P) | in | Counter. | |
j | integer(kind=I4P) | in | Counter. | |
k | integer(kind=I4P) | in | Counter. | |
ngc | integer(kind=I4P) | in | Ghost cells number. | |
nv_c | integer(kind=I4P) | in | Number of conservative varibales in q vector. | |
weno_s | integer(kind=I4P) | in | Weno stencils number/dimension. | |
evmax | real(kind=R8P) | in | Maximum eigenvalue. | |
elw | real(kind=R8P) | in | Left eigenvectors for WENO reconstruction. | |
q | real(kind=R8P) | in | Auxiliary variables. | |
chi | real(kind=R8P) | in | Speed coefficient for D & B div-cleaning. | |
fmpc | real(kind=R8P) | inout | Fluxes -+ 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
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
self | class(prism_cpu_object) | inout | The equation. | |
ivar | integer(kind=I4P) | in | Variable (start) index in q. | |
q | real(kind=R8P) | inout | Field variables. |
Call graph
write_current_behavior_tab
fortran
subroutine write_current_behavior_tab(filename, current_density, time)Arguments
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
filename | character(len=*) | in | ||
current_density | real(kind=R8P) | in | ||
time | real(kind=R8P) | in |
Call graph
write_single_particle_output
fortran
subroutine write_single_particle_output(filename, time, q_pic)Arguments
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
filename | character(len=*) | in | ||
time | real(kind=R8P) | in | ||
q_pic | real(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
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
self | class(prism_cpu_object) | in | The object | |
q | real(kind=R8P) | in | Field variables. | |
n_x | integer(kind=I4P) | in | Number of cells in each direction and number | |
n_y | integer(kind=I4P) | in | Number of cells in each direction and number | |
n_z | integer(kind=I4P) | in | Number of cells in each direction and number | |
n_b | integer(kind=I4P) | in | Number of cells in each direction and number | |
mean_value | real(kind=R8P) | out | Mean 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
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
self | class(prism_cpu_object) | in | The realm. |
Call graph
dotproduct
Compute the scalar (dot) product.
Returns: real(kind=R8P)
fortran
function dotproduct(a, b) result(dot)Arguments
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
a | real(kind=R8P) | in | Left hand side. | |
b | real(kind=R8P) | in | Left hand side. |
Call graph
crossproduct
Returns: real(kind=R8P)
fortran
function crossproduct(a, b) result(cross)Arguments
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
a | real(kind=R8P) | in | Left hand side. | |
b | real(kind=R8P) | in | Left hand side. |
Call graph