Appearance
adam_patch_cpu_object
ADAM, PATCH Poisson solver with Adpative mesh Refinement for HPC computing class definition, CPU backend.
Source: src/app/patch/cpu/adam_patch_cpu_object.F90
Dependencies
Contents
- patch_cpu_object
- allocate_cpu
- finalize
- initialize
- amr_update
- compute_phi
- mark_by_geo
- mark_by_grad_var
- move_phi
- refine_uniform
- integrate_eikonal
- load_restart_files
- save_xh5f
- save_residuals
- save_restart_files
- save_simulation_data
- set_boundary_conditions
- set_initial_conditions
- update_ghost
- integrate
- simulate
Derived Types
patch_cpu_object
Maxwell equations system class definition, CPU backend.
Inheritance
Extends: patch_common_object
Components
| Name | Type | Attributes | Description |
|---|---|---|---|
mpih | type(mpih_object) | MPI handler. | |
adam | type(adam_object) | ADAM. | |
field | type(field_object) | pointer | The field. |
grid | type(grid_object) | pointer | The grid. |
amr | type(amr_object) | AMR marker handler. | |
ib | type(ib_object) | Immersed Boundary (IB) handler. | |
flail | type(flail_object) | Linear algebra methods handler. | |
io | type(patch_io_object) | IO handler. | |
ic | type(patch_ic_object) | Initial Conditions (IC) handler. | |
bc | type(patch_bc_object) | Boundary Conditions (BC) handler. | |
time | type(patch_time_object) | Time handler. | |
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 conservative/primitive variables. |
q | real(kind=R8P) | allocatable | Potential field cell centered variable. |
r | real(kind=R8P) | allocatable | Rho function. |
q_name | character(len=3) | Potential field name. | |
dq | real(kind=R8P) | allocatable | Residuals right hand side. |
Type-Bound Procedures
| Name | Attributes | Description |
|---|---|---|
allocate_common | pass(self) | Allocate common data. |
initialize_common | pass(self) | Initialize the equation common data. |
allocate_cpu | pass(self) | Allocate CPU data. |
finalize | pass(self) | Finalize simulation. |
initialize | pass(self) | Initialize the equation. |
amr_update | pass(self) | Do AMR update. |
compute_phi | pass(self) | Compute phi, distance from IB solid. |
mark_by_geo | pass(self) | Mark blocks to be refined/derefined by a geometric constrain. |
mark_by_grad_var | pass(self) | Mark blocks to be refined/derefined by a grad(var) value. |
move_phi | pass(self) | Move phi. |
refine_uniform | pass(self) | Refine all blocks uniformly. |
integrate_eikonal | pass(self) | Integrate eikonal equation. |
load_restart_files | pass(self) | Load restart files. |
save_xh5f | pass(self) | Save simulation data in XH5F format. |
save_residuals | pass(self) | Save residuals history. |
save_restart_files | pass(self) | Save restart files. |
save_simulation_data | pass(self) | Save all simulation data. |
set_boundary_conditions | pass(self) | Set boundary conditions of equation. |
set_initial_conditions | pass(self) | Set initial conditions (and coils) of equation. |
update_ghost | pass(self) | Update ghost cells and set boundary conditions. |
integrate | pass(self) | Integrate the equation. |
simulate | pass(self) | Perform the simulation. |
Subroutines
allocate_cpu
Allocate CPU data.
fortran
subroutine allocate_cpu(self)Arguments
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
self | class(patch_cpu_object) | inout | The equation. |
Call graph
finalize
Finalize simulation.
fortran
subroutine finalize(self)Arguments
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
self | class(patch_cpu_object) | inout | The equation. |
Call graph
initialize
Initialize the equation.
fortran
subroutine initialize(self, filename)Arguments
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
self | class(patch_cpu_object) | inout | The equation. | |
filename | character(len=*) | in | Input file name. |
Call graph
amr_update
Do AMR update.
fortran
subroutine amr_update(self)Arguments
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
self | class(patch_cpu_object) | inout | The equation. |
Call graph
compute_phi
Compute phi, distance from IB solid.
fortran
subroutine compute_phi(self)Arguments
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
self | class(patch_cpu_object) | inout | The equation. |
Call graph
mark_by_geo
Mark blocks to be refined/derefined by a geometric constrain.
fortran
subroutine mark_by_geo(self, delta_fine, delta_coarse, threshold, do_init)Arguments
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
self | class(patch_cpu_object) | inout | The equation. | |
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 |
Call graph
mark_by_grad_var
Mark blocks to be refined/derefined by a grad(var) value.
fortran
subroutine mark_by_grad_var(self, grad_tol, delta_type, delta_fine, delta_coarse, ivar, threshold, do_init)Arguments
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
self | class(patch_cpu_object) | inout | The equation. | |
grad_tol | real(kind=R8P) | in | Gradiend 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. | |
ivar | integer(kind=I4P) | in | optional | Variable for marking. |
threshold | real(kind=R8P) | in | optional | Threshold for sphere proximity. |
do_init | logical | in | optional | Re-initialize refinements queries. |
Call graph
move_phi
Move phi and the actual tree representation.
fortran
subroutine move_phi(self, velocity, s)Arguments
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
self | class(patch_cpu_object) | inout | The equation. | |
velocity | real(kind=R8P) | in | Velocity of the movement. | |
s | integer(kind=I4P) | in | Solid index. |
Call graph
refine_uniform
Refine all blocks uniformly.
fortran
subroutine refine_uniform(self, refinement_levels)Arguments
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
self | class(patch_cpu_object) | inout | The equation. | |
refinement_levels | integer(kind=I4P) | in | Number of refinement to be performed. |
Call graph
integrate_eikonal
Integrate eikonal equation.
fortran
subroutine integrate_eikonal(self, q)Arguments
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
self | class(patch_cpu_object) | inout | The equation. | |
q | real(kind=R8P) | inout | Conservative variables. |
Call graph
load_restart_files
Save restart files.
fortran
subroutine load_restart_files(self, t, time)Arguments
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
self | class(patch_cpu_object) | inout | The equation. | |
t | integer(kind=I4P) | out | Time iteration. | |
time | real(kind=R8P) | out | Time. |
Call graph
save_xh5f
Save simulation data in HDF5 format.
fortran
subroutine save_xh5f(self, output_basename, with_ghost)Arguments
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
self | class(patch_cpu_object) | inout | The equation. | |
output_basename | character(len=*) | in | optional | Output basename. |
with_ghost | logical | in | optional | Flag to save ghost cells. |
Call graph
save_residuals
Save residuals history.
fortran
subroutine save_residuals(self)Arguments
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
self | class(patch_cpu_object) | inout | The equation. |
Call graph
save_restart_files
Save restart files.
fortran
subroutine save_restart_files(self)Arguments
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
self | class(patch_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(patch_cpu_object) | inout | The equation. |
Call graph
set_boundary_conditions
Set boundary conditions of equation.
fortran
subroutine set_boundary_conditions(self, q)Arguments
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
self | class(patch_cpu_object) | in | The equation. | |
q | real(kind=R8P) | inout | Conservative variables. |
Call graph
set_initial_conditions
Set initial conditions and coils on field.
fortran
subroutine set_initial_conditions(self)Arguments
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
self | class(patch_cpu_object) | inout | The equation. |
Call graph
update_ghost
Update ghost cells. If not specified all steps are perfermod, syncronous computation
fortran
subroutine update_ghost(self, q, step)Arguments
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
self | class(patch_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. |
Call graph
integrate
Integrate the equation.
fortran
subroutine integrate(self, compute_smoothing)Arguments
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
self | class(patch_cpu_object) | inout | The equation. | |
compute_smoothing | procedure(compute_smoothing_interface) | Smoothing procedure. |
Call graph
simulate
Perform the simulation.
fortran
subroutine simulate(self, filename)Arguments
| Name | Type | Intent | Attributes | Description |
|---|---|---|---|---|
self | class(patch_cpu_object) | inout | The equation. | |
filename | character(len=*) | in | Input file name. |
Call graph