Octopus
mix_tests.F90
Go to the documentation of this file.
1!! Copyright (C) 2024 N. Tancogne-Dejean
2!!
3!! This program is free software; you can redistribute it and/or modify
4!! it under the terms of the GNU General Public License as published by
5!! the Free Software Foundation; either version 2, or (at your option)
6!! any later version.
7!!
8!! This program is distributed in the hope that it will be useful,
9!! but WITHOUT ANY WARRANTY; without even the implied warranty of
10!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
11!! GNU General Public License for more details.
12!!
13!! You should have received a copy of the GNU General Public License
14!! along with this program; if not, write to the Free Software
15!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
16!! 02110-1301, USA.
17!!
18
19#include "global.h"
20
22module mix_tests_oct_m
23 use debug_oct_m
25 use global_oct_m
26 use, intrinsic :: iso_fortran_env
28 use mesh_oct_m
29 use mix_oct_m
32 use space_oct_m
33 use types_oct_m
34
35 implicit none
36
37 private
38 public :: &
40
41contains
42
43 subroutine mix_tests_run()
44 integer :: iter, max_iter
45 real(real64), parameter :: tol = 1e-10_real64
46 real(real64), parameter :: solution = 2.0134438978154419300_real64
47 real(real64) :: rhoin(1,1), rhoout(1,1)
48 type(mix_t) :: smix
49 type(mixfield_t), pointer :: mixfield
50 type(derivatives_t) :: der
51 type(mesh_t), target :: mesh
52 type(space_t) :: space
53
54 max_iter = 1000
55
56 ! Fake initializations due to current code design
57 der%dim = 1
58 der%mesh => mesh
59 der%mesh%use_curvilinear = .false.
60 der%mesh%volume_element = m_one
61 der%mesh%parallel_in_domains = .false.
63
64 ! Initialization
65 call mix_init(smix, global_namespace, space, der, 1, 1, func_type_ = type_float)
66 call mix_get_field(smix, mixfield)
67
68 !Starting point
69 rhoin = m_one
70
71 write(message(1),'(a, f20.13)') "Starting point : ", rhoin
72 call messages_info(1)
73
74 do iter = 1, max_iter
75 call mixfield_set_vin(mixfield, rhoin)
76 call test_function(rhoin(1,1), rhoout(1,1))
77
78 write(message(1),'(a, i3, a, f20.13)') "At iter. ", iter, " the value is ", rhoout
79 call messages_info(1)
80
81 ! Check convergence
82 if(abs(rhoout(1,1)-rhoin(1,1)) < tol) exit
83
84 call mixfield_set_vout(mixfield, rhoout)
85 call mixing(global_namespace, smix)
86 call mixfield_get_vnew(mixfield, rhoin)
87 call mixfield_set_vin(mixfield, rhoin)
88 end do
89
90 message(1) = ''
91 if (iter < max_iter) then
92 write(message(2), '(a,f25.15,a,i3,a)') "Computed fixed point ", rhoout, " after ", iter, " iterations."
93 else
94 write(message(2), '(a,f25.15,a,i3,a)') "Failure: fixed point ", rhoout, " after ", iter, " iterations."
95 end if
96 write(message(3),'(a, f25.15, a )') "Numerical fixed point should be ", solution, " ."
97 call messages_info(3)
98
99 end subroutine mix_tests_run
100
104 subroutine test_function(x, fx)
105 real(real64), intent(in) :: x
106 real(real64), intent(out) :: fx
107
108 fx = sin(x) + atan(x)
109
110 end subroutine test_function
111
112end module mix_tests_oct_m
113
114!! Local Variables:
115!! mode: f90
116!! coding: utf-8
117!! End:
double sin(double __x) __attribute__((__nothrow__
double atan(double __x) __attribute__((__nothrow__
This module calculates the derivatives (gradients, Laplacians, etc.) of a function.
real(real64), parameter, public m_one
Definition: global.F90:188
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
subroutine, public messages_info(no_lines, iunit, verbose_limit, stress, all_nodes, namespace)
Definition: messages.F90:624
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:160
subroutine, public mixing(namespace, smix)
Main entry-point to SCF mixer.
Definition: mix.F90:800
subroutine, public mix_get_field(this, mixfield)
Definition: mix.F90:792
subroutine, public mix_init(smix, namespace, space, der, d1, d2, def_, func_type_, prefix_)
Initialise mix_t instance.
Definition: mix.F90:265
This module implements unit tests for the mixing methods.
Definition: mix_tests.F90:115
subroutine, public mix_tests_run()
Definition: mix_tests.F90:137
subroutine test_function(x, fx)
A simple test function.
Definition: mix_tests.F90:198
type(namespace_t), public global_namespace
Definition: namespace.F90:132
type(type_t), public type_float
Definition: types.F90:133