Octopus
convergence_criterion.F90
Go to the documentation of this file.
1!! Copyright (C) 2020 N. Tancogne-Dejean
2!!
3!! This Source Code Form is subject to the terms of the Mozilla Public
4!! License, v. 2.0. If a copy of the MPL was not distributed with this
5!! file, You can obtain one at https://mozilla.org/MPL/2.0/.
6!!
7
8#include "global.h"
9
11 use debug_oct_m
12 use global_oct_m
15 use unit_oct_m
16
17 implicit none
18
19 private
20 public :: &
25
26
27 type, abstract :: convergence_criterion_t
28 private
29 real(real64), public :: tol_abs
30 real(real64), public :: tol_rel
31 character(len=:), allocatable, public :: label
32 type(unit_t), pointer, public :: unit => null()
33
34 real(real64), public :: val_abs
35 real(real64), public :: val_rel
36 real(real64), pointer :: value_diff
37 real(real64), pointer :: norm
38
39 contains
40 procedure :: is_converged => convergence_criterion_is_converged
41 procedure :: set_pointers => criterion_set_quantity_pointers
42 procedure :: write_info => criterion_write_info
44
46 type, extends(linked_list_t) :: criterion_list_t
47 private
48 contains
49 procedure :: add => criterion_list_add_node
50 end type criterion_list_t
51
53 private
54 contains
55 procedure :: get_next => criterion_iterator_get_next
57
58
59contains
60
61 ! ---------------------------------------------------------
62 subroutine criterion_write_info(this, iunit)
63 class(convergence_criterion_t), intent(inout) :: this
64 integer, intent(in) :: iunit
65
66 push_sub(criterion_write_info)
67
68 if (associated(this%unit)) then
69 write(iunit, '(6x,a,a,a,es15.8,a,es15.8,4a)') 'abs_', this%label, ' = ', &
70 units_from_atomic(this%unit, this%val_abs), &
71 ' (', units_from_atomic(this%unit, this%tol_abs), ')', &
72 ' [', trim(units_abbrev(this%unit)), ']'
73 else
74 write(iunit, '(6x,a,a,a,es15.8,a,es15.8,a)') 'abs_', this%label, ' = ', this%val_abs, ' (', this%tol_abs, ')'
75 end if
76 write(iunit, '(6x,a,a,a,es15.8,a,es15.8,a)') 'rel_', this%label, ' = ', this%val_rel, ' (', this%tol_rel, ')'
77
79 end subroutine criterion_write_info
80
81 ! ---------------------------------------------------------
83 subroutine convergence_criterion_is_converged(this, is_converged)
84 class(convergence_criterion_t), intent(inout) :: this
85 logical, intent(out) :: is_converged
86
88
89 this%val_abs = abs(this%value_diff)
90
91 assert(associated(this%norm))
92 if (abs(this%norm) < m_tiny) then
93 this%val_rel = m_huge
94 else
95 this%val_rel = this%val_abs / abs(this%norm)
96 end if
97
98 if (this%tol_abs <= m_zero) then
99 is_converged = .true.
100 else
101 is_converged = this%val_abs < this%tol_abs
102 end if
104 if (this%tol_rel > m_zero) then
105 is_converged = is_converged .and. this%val_rel < this%tol_rel
106 end if
107
108
111
112 ! ---------------------------------------------------------
114 subroutine criterion_set_quantity_pointers(this, value_diff, value_norm)
115 class(convergence_criterion_t), intent(inout) :: this
116 real(real64), target, intent(in) :: value_diff
117 real(real64), target, intent(in) :: value_norm
118
121 this%value_diff => value_diff
122 this%norm => value_norm
126
127 ! ---------------------------------------------------------
129 class(convergence_criterion_t), intent(inout) :: this
132
133 nullify(this%value_diff)
134 nullify(this%norm)
135 nullify(this%unit)
136 if (allocated(this%label)) then
137 ! No safe_deallocate here, because it was not allocated with safe_allocate.
138 deallocate(this%label)
139 end if
140
143
144 ! ---------------------------------------------------------
145 subroutine criterion_list_add_node(this, criterion)
146 class(criterion_list_t) :: this
147 class(convergence_criterion_t), target :: criterion
150
151 call this%add_ptr(criterion)
152
154 end subroutine criterion_list_add_node
156 ! ---------------------------------------------------------
157 function criterion_iterator_get_next(this) result(criterion)
158 class(criterion_iterator_t), intent(inout) :: this
159 class(convergence_criterion_t), pointer :: criterion
160
162
163 select type (ptr => this%get_next_ptr())
164 class is (convergence_criterion_t)
165 criterion => ptr
166 class default
167 assert(.false.)
168 end select
169
171 end function criterion_iterator_get_next
172
174
175!! Local Variables:
176!! mode: f90
177!! coding: utf-8
178!! End:
class(convergence_criterion_t) function, pointer criterion_iterator_get_next(this)
subroutine criterion_set_quantity_pointers(this, value_diff, value_norm)
Setting pointers to the in, out and norm values.
subroutine, public convergence_criterion_end(this)
subroutine convergence_criterion_is_converged(this, is_converged)
Is the convergence reached ?
subroutine criterion_list_add_node(this, criterion)
subroutine criterion_write_info(this, iunit)
real(real64), parameter, public m_huge
Definition: global.F90:205
real(real64), parameter, public m_zero
Definition: global.F90:187
real(real64), parameter, public m_tiny
Definition: global.F90:204
This module implements fully polymorphic linked lists, and some specializations thereof.
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:132
character(len=20) pure function, public units_abbrev(this)
Definition: unit.F90:223
These classes extend the list and list iterator to make a criterion list.
This class implements an iterator for the polymorphic linked list.
This class implements a linked list of unlimited polymorphic values.
int true(void)