All Classes Files Functions Variables Pages
module_lagrange.f90
Go to the documentation of this file.
1 !=====================================================================================================================================
2 ! EFISPEC3D !
3 ! (Elements FInis SPECtraux 3D) !
4 ! !
5 ! http://efispec.free.fr !
6 ! !
7 ! !
8 ! This file is part of EFISPEC3D !
9 ! Please refer to http://efispec.free.fr if you use it or part of it !
10 ! !
11 ! !
12 !1 ---> French License: CeCILL V2 !
13 ! !
14 ! Copyright BRGM 2009 contributeurs : Florent DE MARTIN !
15 ! David MICHEA !
16 ! Philippe THIERRY !
17 ! !
18 ! Contact: f.demartin at brgm.fr !
19 ! !
20 ! Ce logiciel est un programme informatique servant a resoudre l'equation du !
21 ! mouvement en trois dimensions via une methode des elements finis spectraux. !
22 ! !
23 ! Ce logiciel est regi par la licence CeCILL soumise au droit francais et !
24 ! respectant les principes de diffusion des logiciels libres. Vous pouvez !
25 ! utiliser, modifier et/ou redistribuer ce programme sous les conditions de la !
26 ! licence CeCILL telle que diffusee par le CEA, le CNRS et l'INRIA sur le site !
27 ! "http://www.cecill.info". !
28 ! !
29 ! En contrepartie de l'accessibilite au code source et des droits de copie, de !
30 ! modification et de redistribution accordes par cette licence, il n'est offert !
31 ! aux utilisateurs qu'une garantie limitee. Pour les memes raisons, seule une !
32 ! responsabilite restreinte pese sur l'auteur du programme, le titulaire des !
33 ! droits patrimoniaux et les concedants successifs. !
34 ! !
35 ! A cet egard l'attention de l'utilisateur est attiree sur les risques associes !
36 ! au chargement, a l'utilisation, a la modification et/ou au developpement et a !
37 ! la reproduction du logiciel par l'utilisateur etant donne sa specificite de !
38 ! logiciel libre, qui peut le rendre complexe a manipuler et qui le reserve donc !
39 ! a des developpeurs et des professionnels avertis possedant des connaissances !
40 ! informatiques approfondies. Les utilisateurs sont donc invites a charger et !
41 ! tester l'adequation du logiciel a leurs besoins dans des conditions permettant !
42 ! d'assurer la securite de leurs systemes et ou de leurs donnees et, plus !
43 ! generalement, a l'utiliser et l'exploiter dans les memes conditions de !
44 ! securite. !
45 ! !
46 ! Le fait que vous puissiez acceder a cet en-tete signifie que vous avez pris !
47 ! connaissance de la licence CeCILL et que vous en avez accepte les termes. !
48 ! !
49 ! !
50 !2 ---> International license: GNU GPL V3 !
51 ! !
52 ! EFISPEC3D is a computer program that solves the three-dimensional equations of !
53 ! motion using a finite spectral-element method. !
54 ! !
55 ! Copyright (C) 2009 Florent DE MARTIN !
56 ! !
57 ! Contact: f.demartin at brgm.fr !
58 ! !
59 ! This program is free software: you can redistribute it and/or modify it under !
60 ! the terms of the GNU General Public License as published by the Free Software !
61 ! Foundation, either version 3 of the License, or (at your option) any later !
62 ! version. !
63 ! !
64 ! This program is distributed in the hope that it will be useful, but WITHOUT ANY !
65 ! WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A !
66 ! PARTICULAR PURPOSE. See the GNU General Public License for more details. !
67 ! !
68 ! You should have received a copy of the GNU General Public License along with !
69 ! this program. If not, see http://www.gnu.org/licenses/. !
70 ! !
71 ! !
72 !3 ---> Third party libraries !
73 ! !
74 ! EFISPEC3D uses the following of third party libraries: !
75 ! !
76 ! --> METIS 5.1.0 !
77 ! see http://glaros.dtc.umn.edu/gkhome/metis/metis/overview !
78 ! !
79 ! --> Lib_VTK_IO !
80 ! see S. Zaghi's website: https://github.com/szaghi/Lib_VTK_IO !
81 ! !
82 ! --> INTERP_LINEAR !
83 ! see J. Burkardt website: http://people.sc.fsu.edu/~jburkardt/ !
84 ! !
85 ! --> SAC !
86 ! http://ds.iris.edu/files/sac-manual/ !
87 ! !
88 ! --> EXODUS II !
89 ! http://sourceforge.net/projects/exodusii/ !
90 ! !
91 ! --> NETCDF !
92 ! http://www.unidata.ucar.edu/software/netcdf/ !
93 ! !
94 ! --> HDF5 !
95 ! http://www.hdfgroup.org/HDF5/ !
96 ! !
97 ! Some of these libraries are located in directory lib and pre-compiled !
98 ! with intel compiler for x86_64 platform. !
99 ! !
100 ! !
101 !4 ---> Related Articles !
102 ! !
103 ! De Martin, F., Matsushima, M., Kawase, H. (BSSA, 2013) !
104 ! Impact of geometric effects on near-surface Green's functions !
105 ! doi:10.1785/0120130039 !
106 ! !
107 ! Aochi, H., Ducellier, A., Dupros, F., Delatre, M., Ulrich, T., De Martin, F., !
108 ! Yoshimi, M., (Pure Appl. Geophys. 2013) !
109 ! Finite Difference Simulations of Seismic Wave Propagation for the 2007 Mw 6.6 !
110 ! Niigata-ken Chuetsu-Oki Earthquake: Validity of Models and Reliable !
111 ! Input Ground Motion in the Near-Field !
112 ! doi:10.1007/s00024-011-0429-5 !
113 ! !
114 ! De Martin, F. (BSSA, 2011) !
115 ! Verification of a Spectral-Element Method Code for the Southern California !
116 ! Earthquake Center LOH.3 Viscoelastic Case !
117 ! doi:10.1785/0120100305 !
118 ! !
119 !=====================================================================================================================================
120 
123 
127 
128  use mpi
129 
130  implicit none
131 
132  private
133 
134  public :: lagrap
135  public :: lagrap_geom
136  public :: lagrad
137  public :: lagrad_geom
138  public :: hexa_lagrange_interp
139  public :: quad_lagrange_interp
140 
141  contains
142 
143 !
144 !
149 !***********************************************************************************************************************************************************************************
150  real function lagrap(i,x,n)
151 !***********************************************************************************************************************************************************************************
152 
153  use mod_global_variables, only :&
154  rg_gll_abscissa&
155  ,rg_gll_abscissa_dist
156 
157  implicit none
158 
159  real , intent(in) :: x
160  integer, intent(in) :: i
161  integer, intent(in) :: n
162 
163  integer :: j
164 
165  lagrap = 1.0
166  do j = 1,n
167  if (j.ne.i) lagrap = lagrap*( (x-rg_gll_abscissa(j))*rg_gll_abscissa_dist(j,i) )
168  enddo
169 
170 !***********************************************************************************************************************************************************************************
171  end function lagrap
172 !***********************************************************************************************************************************************************************************
173 
174 !
175 !
180 !***********************************************************************************************************************************************************************************
181  real function lagrap_geom(i,x,n)
182 !***********************************************************************************************************************************************************************************
183 
184  use mod_global_variables, only :&
185  rg_gnode_abscissa_dist&
186  ,rg_gnode_abscissa
187 
188  implicit none
189 
190  real , intent(in) :: x
191  integer, intent(in) :: i
192  integer, intent(in) :: n
193 
194  integer :: j
195 
196  lagrap_geom = 1.0
197  do j = 1,n
198  if (j.ne.i) lagrap_geom = lagrap_geom*( (x-rg_gnode_abscissa(j))*rg_gnode_abscissa_dist(j,i) )
199  enddo
200 
201 !***********************************************************************************************************************************************************************************
202  end function lagrap_geom
203 !***********************************************************************************************************************************************************************************
204 
205 !
206 !
212 !***********************************************************************************************************************************************************************************
213  real function lagrad(i,x,n)
214 !***********************************************************************************************************************************************************************************
215 
216  use mod_global_variables, only : rg_gll_abscissa&
217  ,rg_gll_abscissa_dist
218 
219  implicit none
220 
221  real , intent(in) :: x
222  integer, intent(in) :: i
223  integer, intent(in) :: n
224 
225  real :: ltmp
226  real :: lprod
227  integer :: j
228  integer :: k
229 
230  lagrad = 0.0
231  lprod = 0.0
232  do j = 1,n
233  if (j.ne.i) then
234  lprod = 1.0
235  do k = 1,n
236  if (k.ne.i) then
237  if(k.eq.j) then
238  ltmp = rg_gll_abscissa_dist(j,i)
239  else
240  ltmp = (x-rg_gll_abscissa(k))*rg_gll_abscissa_dist(k,i)
241  endif
242  lprod = lprod*ltmp
243  endif
244  enddo
245  endif
246  lagrad = lagrad + lprod
247  lprod = 0.0
248  enddo
249 !
250 !***********************************************************************************************************************************************************************************
251  end function lagrad
252 !***********************************************************************************************************************************************************************************
253 
254 !
255 !
261 !***********************************************************************************************************************************************************************************
262  real function lagrad_geom(i,x,n)
263 !***********************************************************************************************************************************************************************************
264 
265  use mod_global_variables, only :&
266  rg_gnode_abscissa&
267  ,rg_gnode_abscissa_dist
268 
269  implicit none
270 
271  real , intent(in) :: x
272  integer, intent(in) :: i
273  integer, intent(in) :: n
274 
275  real :: ltmp
276  real :: lprod
277  integer :: j
278  integer :: k
279 
280  lagrad_geom = 0.0
281  lprod = 0.0
282  do j = 1,n
283  if (j.ne.i) then
284  lprod = 1.0
285  do k = 1,n
286  if (k.ne.i) then
287  if(k.eq.j) then
288  ltmp = rg_gnode_abscissa_dist(j,i)
289  else
290  ltmp = (x-rg_gnode_abscissa(k))*rg_gnode_abscissa_dist(k,i)
291  endif
292  lprod = lprod*ltmp
293  endif
294  enddo
295  endif
296  lagrad_geom = lagrad_geom + lprod
297  lprod = 0.0
298  enddo
299 !
300 !***********************************************************************************************************************************************************************************
301  end function lagrad_geom
302 !***********************************************************************************************************************************************************************************
303 
304 !
305 !
313 !***********************************************************************************************************************************************************************************
314  subroutine hexa_lagrange_interp(gll_values,lag,x,y,z)
315 !***********************************************************************************************************************************************************************************
316 
317  use mod_global_variables, only :&
318  ig_ngll&
319  ,ig_ndof
320 
321  implicit none
322 
323  real, intent( in), dimension(IG_NDOF,IG_NGLL,IG_NGLL,IG_NGLL) :: gll_values
324  real, intent( in), dimension( IG_NGLL,IG_NGLL,IG_NGLL) :: lag
325  real, intent(out) :: x
326  real, intent(out) :: y
327  real, intent(out) :: z
328 
329  integer :: k
330  integer :: l
331  integer :: m
332 
333 
334  x = 0.0
335  y = 0.0
336  z = 0.0
337 
338  do k = 1,ig_ngll
339  do l = 1,ig_ngll
340  do m = 1,ig_ngll
341 
342  x = x + gll_values(1,m,l,k)*lag(m,l,k)
343  y = y + gll_values(2,m,l,k)*lag(m,l,k)
344  z = z + gll_values(3,m,l,k)*lag(m,l,k)
345 
346  enddo
347  enddo
348  enddo
349 
350  return
351 !***********************************************************************************************************************************************************************************
352  end subroutine hexa_lagrange_interp
353 !***********************************************************************************************************************************************************************************
354 
355 !
356 !
364 !***********************************************************************************************************************************************************************************
365  subroutine quad_lagrange_interp(gll_values,lag,x,y,z)
366 !***********************************************************************************************************************************************************************************
367 
368  use mod_global_variables, only :&
369  ig_ngll&
370  ,ig_ndof
371 
372  implicit none
373 
374  real, intent( in), dimension(IG_NDOF,IG_NGLL,IG_NGLL) :: gll_values
375  real, intent( in), dimension( IG_NGLL,IG_NGLL) :: lag
376  real, intent(out) :: x
377  real, intent(out) :: y
378  real, intent(out) :: z
379 
380  integer :: k
381  integer :: l
382 
383 
384  x = 0.0
385  y = 0.0
386  z = 0.0
387 
388  do k = 1,ig_ngll
389  do l = 1,ig_ngll
390 
391  x = x + gll_values(1,l,k)*lag(l,k)
392  y = y + gll_values(2,l,k)*lag(l,k)
393  z = z + gll_values(3,l,k)*lag(l,k)
394 
395  enddo
396  enddo
397 
398  return
399 !***********************************************************************************************************************************************************************************
400  end subroutine quad_lagrange_interp
401 !***********************************************************************************************************************************************************************************
402 
403 end module mod_lagrange
real function, public lagrap(i, x, n)
function to compute value of order Lagrange polynomial of the GLL node i at abscissa : ...
real function, public lagrad_geom(i, x, n)
function to compute value of the derivative of order Lagrange polynomial of geometric node i at absc...
This module defines all global variables of EFISPEC3D. Scalar variables are initialized directly in t...
real function, public lagrap_geom(i, x, n)
function to compute value of order Lagrange polynomial of geometric node i at abscissa : ...
subroutine, public quad_lagrange_interp(gll_values, lag, x, y, z)
This subroutine interpolates GLL nodes -values of a quadrangle element using Lagrange polynomials...
This module contains functions to compute Lagrange polynomials and its derivatives; and to interpolat...
real function, public lagrad(i, x, n)
function to compute value of the derivative of order Lagrange polynomial of GLL node i at abscissa ...
subroutine, public hexa_lagrange_interp(gll_values, lag, x, y, z)
This subroutine interpolates GLL nodes -values of a hexahedron element using Lagrange polynomials...