Usage examples

Read a single frame (indexes.f90)

program indexes_
    use iso_fortran_env, only: real64, int64
    use chemfiles
    implicit none
    ! Chemfiles types declaration uses the "chfl_" prefix
    type(chfl_trajectory) :: trajectory
    type(chfl_frame) :: frame

    real(real64), dimension(:, :), pointer :: positions
    integer(int64), dimension(:), allocatable :: indexes
    integer(int64) :: natoms, i, j
    integer :: status

    call frame%init()
    call trajectory%open("filename.xyz", "r")

    call trajectory%read(frame, status=status)
    if (status /= 0) then
        stop "Error"
    end if

    call frame%positions(positions, natoms)
    allocate(indexes(natoms))

    indexes = -1
    j = 1
    do i=1,natoms
        if (positions(1, i) < 5) then
            indexes(j) = i
            j = j + 1
        end if
    end do

    write(*,*) "Atoms with x < 5: "
    do i=1,natoms
        if (indexes(i) == -1) then
            exit
        end if
        write(*,*) "  - ", indexes(i)
    end do

    ! Cleanup the allocated memory
    deallocate(indexes, positions)
    call trajectory%close()
    call frame%free()
end program

Read Multiple frames (rmsd.f90)

program rmsd_
    use iso_fortran_env, only: real64, int64
    use chemfiles
    implicit none

    type(chfl_trajectory) :: trajectory
    type(chfl_frame)      :: frame
    real(real64), dimension(:, :), pointer :: positions
    real(real64), dimension(:),    allocatable :: distances
    integer(int64) :: nsteps = 0, natoms=0, i
    real(real64) :: distance = 0, mean = 0, rmsd = 0
    integer :: status

    call trajectory%open("filename.nc", "r", status=status)
    if (status /= 0) stop "Error while opening input file"
    call frame%init()

    call trajectory%nsteps(nsteps)
    allocate(distances(nsteps))

    ! Accumulate the distances to the origin of the 10th atom throughtout the
    ! trajectory
    do i=1,nsteps
        call trajectory%read(frame, status)
        if (status /= 0) stop "Error while reading the frame"

        ! Only allocate on the first iteration. This assume a constant number
        ! of particles
        if (i == 0) then
            call frame%atoms_count(natoms)
            allocate(positions(3, natoms))
        end if

        ! Position of the 10th atom
        call frame%positions(positions, natoms, status)
        if (status /= 0) stop "Error while getting the positions"
        distance = sqrt(positions(1, 10)*positions(1, 10) + &
                        positions(2, 10)*positions(2, 10) + &
                        positions(3, 10)*positions(3, 10))
        distances(i) = distance
    end do

    do i=1,nsteps
        mean = mean + distances(i)
    end do
    mean = mean / nsteps

    do i=1,nsteps
        rmsd = rmsd + (mean - distances(i))*(mean - distances(i));
    end do
    rmsd = rmsd / nsteps;
    rmsd = sqrt(rmsd);

    write(*, *) "Root-mean square displacement is: ", rmsd

    ! Cleanup the allocated memory
    call trajectory%close()
    call frame%free()
    deallocate(distances, positions)
end program

Write frames (convert.f90)

program convert
    use iso_fortran_env, only: int64
    use chemfiles
    implicit none

    type(chfl_trajectory) :: input, ouput_file
    type(chfl_frame) :: frame
    type(chfl_cell) :: cell
    type(chfl_topology) :: water_topology
    type(chfl_atom) :: O, H

    integer(int64) :: nsteps, i
    integer :: status

    call input%open("water.xyz", "r", status=status)
    if (status /= 0) stop "Error while opening input file"
    call ouput_file%open("water.pdb", "w")

    call frame%init()
    call water_topology%init()
    ! Orthorombic UnitCell with lengths of 20, 15 and 35 A
    call cell%init([20d0, 15d0, 35d0])

    ! Create Atoms
    call O%init("O")
    call H%init("H")

    ! Fill the topology with one water molecule
    call water_topology%add_atom(O)
    call water_topology%add_atom(H)
    call water_topology%add_atom(H)
    call water_topology%add_bond(0_int64, 1_int64)
    call water_topology%add_bond(0_int64, 2_int64)

    call input%nsteps(nsteps)

    do i=1,nsteps
        call input%read(frame, status)
        if (status /= 0) stop "Error while reading input file"
        ! Set the frame cell and topology
        call frame%set_cell(cell)
        call frame%set_topology(water_topology)
        ! Write the frame to the output file, using PDB format
        call ouput_file%write(frame, status)
        if (status /= 0) stop "Error while writing ouput file"
    end do

    ! Cleanup the allocated memory
    call input%close()
    call ouput_file%close()
    call frame%free()
    call cell%free()
    call water_topology%free()
    call O%free()
    call H%free()
end program