ChiDG APIΒΆ

The ChiDG API is a Fortran Class(chidg_t) that provides an interface to the ChiDG library.

type :: chidg_t

    type(chidg_data_t)                      :: data
    class(time_integrator_t),   allocatable :: time_integrator
    class(nonlinear_solver_t),  allocatable :: nonlinear_solver
    class(linear_solver_t),     allocatable :: linear_solver
    class(preconditioner_t),    allocatable :: preconditioner

contains

    ! Open/Close
    procedure   :: start_up
    procedure   :: shut_down

    ! Options
    procedure   :: set

    ! IO
    procedure   :: read_mesh
    procedure   :: read_fields
    procedure   :: write_mesh
    procedure   :: write_fields

    ! Run
    procedure   :: run

end type chidg_t

Example:

! This is an example of how a driver might use the ChiDG API to run an
! analysis. This is similar to the file src/driver.f90 that is used
! to create the ChiDG executable.

! Import the ChiDG API class
use type_chidg, only: chidg_t

! Create an instance of the ChiDG API
type(chidg_t)   :: chidg

! Initialize MPI and ChiDG infrastructure
call chidg%start_up('mpi')
call chidg%start_up('core')

! Set algorithms
call chidg%set('Time Integrator' , algorithm='Steady'      )
call chidg%set('Nonlinear Solver', algorithm='Quasi-Newton')
call chidg%set('Linear Solver'   , algorithm='FGMRES'      )
call chidg%set('Preconditioner'  , algorithm='ILU0'        )
call chidg%set('Solution Order'  , integer_input=2         )

! Read grid (solution if available)
call chidg%read_mesh('grid.h5')
call chidg%read_fields('grid.h5')

! Run ChiDG analysis
call chidg%run()

API procedures:

chidg%start_up('option')

Performs start-up routines on a ChiDG instance. Options are:

mpi
Calls MPI_Init, initializes IRANK, NRANK
core
Register functions, equations, models, operators, bcs, compute reference
matrices. Initialize ChiDG_COMM MPI communicator. Default communicator
is MPI_COMM_WORLD.
namelist
Read namelist input file, chidg.nml
chidg%shut_down('option')

Performs shut-down routines on a ChiDG instance. Options are:

mpi
Calls MPI_Finalize.
core
Close log files, close HDF interface.
log
Close log files.
chidg%set('String Selector', algorithm='string', integer_index=value)

Set algorithms and parameters for the ChiDG instance.

String Selector Expected input
Time Integrator algorithm=’string’
Nonlinear Solver algorithm=’string’
Linear Solver algorithm=’string’
Preconditioner algorithm=’string’
Solution Order integer_input=value
chidg%read_mesh('file name')

Read domains and boundary conditions from the ChiDG-formatted HDF file indicated by the string parameter passed in.

Note

All calls to chidg%set() shall occur before this.

chidg%read_fields('file name')

Read solution from the ChiDG-formatted HDF file indicated by the string parameter passed in.

Note

A grid shall already have been read/initialized before a call to chidg%read_fields('file').

chidg%write_mesh('file name')

Write domains and boundary conditions to the ChiDG-formatted HDF file indicated by the string parameter passed in. It ‘file name’ does not exist, a new file will be created.

chidg%write_fields('file name')

Write primary fields to the ChiDG-formatted HDF file indicated by the string parameter passed in.

chidg%run()

Run the ChiDG analysis. This calls a certain number of ‘steps’ on the time integrator.

Note

All reading/initializing/setting shall occur before this.