Choral
stimulation_mod.f90
Go to the documentation of this file.
1 
27 
29 
31  use basic_tools
32  use real_type
34  use funclib
35 
36  implicit none
37 
38  public :: stimulation
39 
40 
41  ! %----------------------------------------%
42  ! | |
43  ! | DERIVED TYPE |
44  ! | |
45  ! %----------------------------------------%
70 
72  procedure(rtor), nopass, pointer :: f => f_c3
73 
75  real(RP), dimension(3) :: x0 = 0._rp
76 
78  real(RP) :: t0 = 3._rp
79 
81  real(RP) :: rt = 1._rp
82 
84  real(RP) :: rx = 0.15_rp
85 
87  real(RP) :: l = 1.0_rp
88 
89  contains
90 
91  procedure :: I_app => stim
92 
93  end type stimulation
94 
95 
96  ! %----------------------------------------%
97  ! | |
98  ! | GENERIc SUBROUTINES |
99  ! | |
100  ! %----------------------------------------%
101  interface stimulation
102  module procedure stimulation_create
103  end interface stimulation
104 
105 
106 contains
107 
123  function stimulation_create(f, L, t0, Rt, Rx, x0) result(s)
124  type(stimulation) :: s
125 
126  procedure(rtor) :: F
127  real(RP), intent(in) :: L
128  real(RP), intent(in) :: t0
129  real(RP), intent(in) :: Rt
130  real(RP), intent(in) :: Rx
131  real(RP), dimension(3), optional, intent(in) :: x0
132 
133  s%t0 = t0
134  s%Rt = rt
135  s%Rx = rx
136  s%L = l
137  s%F => f
138 
139  if (present(x0)) then
140  s%x0 = x0
141  else
142  s%x0 = 0._rp
143  end if
144 
145  if (.NOT.associated(s%f)) call quit("stimulation_mod:&
146  & stimulation_create: invalid argument 'f'")
147  if (.NOT.(s%Rt>0._rp)) call quit("stimulation_mod:&
148  & stimulation_create: invalid argument 'Rt'")
149  if (.NOT.(s%Rx>0._rp)) call quit("stimulation_mod:&
150  & stimulation_create: invalid argument 'Rx'")
151 
152  end function stimulation_create
153 
154 
155  function stim(s, x, t) result(res)
156  class(stimulation) , intent(in) :: s
157  real(RP) :: res
158  real(RP), dimension(3), intent(in) :: x
159  real(RP) , intent(in) :: t
160 
161  real(RP) :: t2, r
162 
163  t2 = abs( t - s%t0)/s%Rt
164 
165  r = sqrt(sum( (x-s%x0)**2 ))/s%Rx
166 
167  res = s%f(t2) * s%f(r) * s%L
168 
169  end function stim
170 
171 
172 end module stimulation_mod
173 
real(rp) function stim(s, x, t)
DEFINITION OF ABSTRACT_INTERFACES FOR THE LIBRARY CHORAL
BASIC TOOLS
Definition: basic_tools.F90:8
real(rp) function, public f_c3(t)
C3 function with support [-1,1].
Definition: funcLib.f90:164
subroutine, public quit(message)
Stop program execution, display an error messahe.
PRE-DEFINED NUMERIC FUNCTIONS
Definition: funcLib.f90:25
REAL NUMBERS PRECISION IN CHORAL: selects simple/double/quad
Definition: real_type.F90:33
Derived type stimulation.
type(stimulation) function stimulation_create(f, L, t0, Rt, Rx, x0)
Constructor for the type stimulation
CHORAL CONSTANTS
real(rp) function f(x)
right hand side 'f'
DERIVED TYPE stimulation