The Doxygen
documentation for the main branch of Octopus can be found here.
It is automatically generated at a push to the main branch.
Currently, the links to the source code are broken and will most likely point to a wrong part of the file.
Running doxygen locally
To build the Doxygen
pages locally, go into the folder doc/doxygen/
inside the Octopus tree, and run the command
./create_documentation.sh
This might take a few minutes to complete. In case you need to rebuild the pages more often while writing new documentation, you might want to disable the creation of call graphs, by hanging the line
HAVE_DOT = YES
to
HAVE_DOT = NO
in the file Doxyfile
.
Writing documentation
There are several ways how to include comments into the Fortran source, which are shown in the following code snippet:
In particular, we would like to have documentation of modules, classes, and subroutines/functions, together with their data members and arguments.
General remarks about Doxygen
Doxygen
markers in Fortran are !>, !< and !!. The first two denote the beginning of a Doxygen
comment,
and the direction of the ‘arrow’ indicates what the comment refers to. In general, if a comment is before the
quantity which should be documented, one should use !>, while if the comment is at the end of the line
(e.g. for an function argument), one should use !<. Continuation lines use the !! marker.
Using the wrong marker can lead to wrong and confusing documentation.
A continuation of a comment to a function argument, written as
integer,intent(in)::dim!< Dimension of space
!! needs to be in the range 1 to 4.
is not accepted by findent-octopus
and will lead to a failing test.
Instead, this should be written as
integer,intent(in)::dim!< Dimension of space
!! needs to be in the range 1 to 4.
Markdown and Equations in Doxygen
Doxygen
supports Markdown syntax for formatting text. More details can be found here.
Furthermore, we enabled MathJax, which allows to use LaTeX formatted equations in the code documentation.
Mathematical formulas have to be enclosed by a pair of \f$ for inline formulas, or by the pair of \f[ and \f]
for formulas to be displayed as separate line. See https://www.doxygen.nl/manual/formulas.html.
Module documentation
Before the declaration of the module, there should be a block
!> @brief This modules implements ...
!!
!! Some more details about the modules, how it is related to the code, general theory,
!! references to papers, etc.
!!
module module1_oct_m...end module
There is one peculiarity about the Doxygen usage in Octopus, which is that we need to pass the sources through
the cpp preprocessor in order to handle some macros, in particular out poor-man’s templating. This has the
side effect that the preprocessor removes all C/C++ like comments from the code, including the ‘//’ in URLs.
URLs can be included in Doxygen comments by using the HTTP and HTTPS macros.
!> @brief The Octopus code
!!
!! The web documentation can be accessed \HTTPS{here,octopus-code.org}.
!!
Class documentation
Before the declaration of the module, there should be a block
!> @brief This class implements ..
!!
!! Some more details about the class, which is not in the module description
!! or in the documentation of member data or procedure
!!
type class1_tprivate
class(system_t)::system!< This system is related to the current object
logical,public::is_initialized!< indicates whether all initialization has been successfully done
...contains! Below are the methods, the class implements
procedure::do_something=>class1_do_something!< @copydoc module1_oct_m::class1_do_something
procedure::finished=>class1_finished!< @copydoc module1_oct_m::class1_finished
...end type
The @copydoc command here instructs
Doxygen
to copy the documentation from the actual function definition for the class method.
Note that the corresponding module name has to be prepended (in C++ codes, this would be the namespace).
Subroutine / Function documentation
The documentation of functions and subroutines is similar to that of classes.
Example
Example
!! Copyright (C) 2019 N. Tancogne-Dejean
!! Copyright (C) 2020 M. Oliveira
!!
!! This program is free software; you can redistribute it and/or modify
!! it under the terms of the GNU General Public License as published by
!! the Free Software Foundation; either version 2, or (at your option)
!! any later version.
!!
!! This program is distributed in the hope that it will be useful,
!! but WITHOUT ANY WARRANTY; without even the implied warranty of
!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
!! GNU General Public License for more details.
!!
!! You should have received a copy of the GNU General Public License
!! along with this program; if not, write to the Free Software
!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
!! 02110-1301, USA.
!!
#include "global.h"
!> @brief This module implements the basic propagator framework.
!!
module propagator_oct_muse algorithm_oct_muse clock_oct_muse debug_oct_muse global_oct_muse messages_oct_muse namespace_oct_muse profiling_oct_muse system_oct_mimplicit none
private
public::&propagator_t!> @brief Abstract class implementing propagators
!!
!! Propagators are implemented as a state machine. This abstract class defines the steps which
!! are independent of the actual propatation algorithm, and also independent of the system, to
!! which is is applied.
type,extends(algorithm_t),abstract::propagator_tprivate
class(system_t),pointer,public::system
type(algorithm_iterator_t)::scf_start!< Options related to predictor-corrector propagators
logical,public::predictor_corrector=.false.integer,public::scf_countinteger,public::max_scf_countFLOAT,public::scf_tolFLOAT,public::final_time=M_ZEROcontains! Below are the list of operations that needs to be implemented
procedure::do_operation=>propagator_do_operation!< @copydoc propagator_oct_m::propagator_do_operation
procedure::finished=>propagator_finished!< @copydoc propagator_oct_m::propagator_finished
procedure::save_scf_start=>propagator_save_scf_start!< @copydoc propagator_oct_m::propagator_save_scf_start
procedure::rewind_scf_loop=>propagator_rewind_scf_loop!< @copydoc propagator_oct_m::propagator_rewind_scf_loop
end type propagator_t!# doc_start general_propagation_operations
! Known propagation operations
character(len=ALGO_LABEL_LEN),public,parameter::&START_SCF_LOOP='START_SCF_LOOP',&END_SCF_LOOP='END_SCF_LOOP',&STORE_CURRENT_STATUS='STORE_CURRENT_STATUS'type(algorithmic_operation_t),public,parameter::&OP_START_SCF_LOOP=algorithmic_operation_t(START_SCF_LOOP,'Starting SCF loop'),&OP_END_SCF_LOOP=algorithmic_operation_t(END_SCF_LOOP,'End of SCF iteration'),&OP_STORE_CURRENT_STATUS=algorithmic_operation_t(STORE_CURRENT_STATUS,'Store current status')!# doc_end
contains!> @brief perform one operation of the state machine
!!
!! this routine performs operations which are general (not system specific).
!!
logical function propagator_do_operation(this,operation)result(done)class(propagator_t),intent(inout)::thistype(algorithmic_operation_t),intent(in)::operationdone=.true.select case(operation%id)case(START_SCF_LOOP)ASSERT(this%predictor_corrector)call this%save_scf_start()this%inside_scf=.true.this%system%accumulated_loop_ticks=0if(debug%info)then
write(message(1),'(a,i3,a)')"Debug: -- SCF iter ",this%scf_count," for '"+trim(this%system%namespace%get())+"'"call messages_info(1,namespace=this%system%namespace)end if
case(END_SCF_LOOP)! Here we first check if we did the maximum number of steps.
! Otherwise, we need check the tolerance
if(this%scf_count==this%max_scf_count)then
if(debug%info)then
message(1)="Debug: -- Max SCF Iter reached for '"+trim(this%system%namespace%get())+"'"call messages_info(1,namespace=this%system%namespace)end if
this%inside_scf=.false.call this%next()else! We reset the pointer to the beginning of the scf loop
if(this%system%is_tolerance_reached(this%scf_tol))then
if(debug%info)then
message(1)="Debug: -- SCF tolerance reached for '"+trim(this%system%namespace%get())+"'"call messages_info(1,namespace=this%system%namespace)end if
this%inside_scf=.false.call this%next()else! We rewind the instruction stack
call this%rewind_scf_loop()! We reset the clocks
call this%system%reset_clocks(this%system%accumulated_loop_ticks)this%system%accumulated_loop_ticks=0if(debug%info)then
write(message(1),'(a,i3,a,a)')"Debug: -- SCF iter ",this%scf_count," for '"+trim(this%system%namespace%get()),"'"call messages_info(1,namespace=this%system%namespace)end if
end if
end if
case defaultdone=.false.end select
end function propagator_do_operation! ---------------------------------------------------------
!> @brief indicate whether a propagation has reached the final time
!!
logical function propagator_finished(this)class(propagator_t),intent(in)::thistype(clock_t)::clock_clock_=this%system%clock+CLOCK_TICKpropagator_finished=clock_%time()>this%final_timeend function propagator_finished!> @brief Save the current iteration state (START_SCF_LOOP) and move to next step
!!
subroutine propagator_save_scf_start(this)class(propagator_t),intent(inout)::thisPUSH_SUB(propagator_save_scf_start)this%scf_start=this%itercall this%next()this%scf_count=0POP_SUB(propagator_save_scf_start)end subroutine propagator_save_scf_start! ---------------------------------------------------------
!> @brief Reset the iteration state to the beginning of the loop (START_SCF_LOOP) and move to next step
!!
subroutine propagator_rewind_scf_loop(this)class(propagator_t),intent(inout)::thisPUSH_SUB(propagator_rewind_scf_loop)this%iter=this%scf_startcall this%next()this%scf_count=this%scf_count+1POP_SUB(propagator_rewind_scf_loop)end subroutine propagator_rewind_scf_loopend module propagator_oct_m!! Local Variables:
!! mode: f90
!! coding: utf-8
!! End: