SPEED
GET_MONITOR_VALUE.f90
Go to the documentation of this file.
1! Copyright (C) 2012 The SPEED FOUNDATION
2! Author: Ilario Mazzieri
3!
4! This file is part of SPEED.
5!
6! SPEED is free software; you can redistribute it and/or modify it
7! under the terms of the GNU Affero General Public License as
8! published by the Free Software Foundation, either version 3 of the
9! License, or (at your option) any later version.
10!
11! SPEED is distributed in the hope that it will be useful, but
12! WITHOUT ANY WARRANTY; without even the implied warranty of
13! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
14! Affero General Public License for more details.
15!
16! You should have received a copy of the GNU Affero General Public License
17! along with SPEED. If not, see <http://www.gnu.org/licenses/>.
18
30
31
32 subroutine get_monitor_value(nb_nod, xq, val, xref, yref, zref, res)
33
34 implicit none
35
36 integer*4 :: nb_nod
37 integer*4 :: i,j,k,h
38
39 real*8 :: xref,yref,zref,res
40 real*8 :: fx,fy,fz
41
42 real*8, dimension(nb_nod) :: xq
43
44 real*8, dimension(nb_nod,nb_nod,nb_nod) :: val
45
46 res = 0.0d0
47
48 do k = 1,nb_nod
49 fz = 1.0d0
50 do h = 1,nb_nod
51 if (h.ne.k) fz = fz * (zref - xq(h)) / (xq(k) - xq(h))
52 enddo
53
54 do j = 1,nb_nod
55 fy = 1.0d0
56 do h = 1,nb_nod
57 if (h.ne.j) fy = fy * (yref - xq(h)) / (xq(j) - xq(h))
58 enddo
59
60 do i = 1,nb_nod
61 fx = 1.0d0
62 do h = 1,nb_nod
63 if (h.ne.i) fx = fx * (xref - xq(h)) / (xq(i) - xq(h))
64 enddo
65
66 res = res + val(i,j,k) * fx * fy * fz
67 enddo
68 enddo
69 enddo
70
71 end subroutine get_monitor_value
72
subroutine get_monitor_value(nb_nod, xq, val, xref, yref, zref, re
Computes the mean value.