My Project
particle.f90
Go to the documentation of this file.
1 
2 
3 
4 
5 
6 
7 
8 
9 
10 
11 
12 !/===========================================================================/
13 ! Copyright (c) 2007, The University of Massachusetts Dartmouth
14 ! Produced at the School of Marine Science & Technology
15 ! Marine Ecosystem Dynamics Modeling group
16 ! All rights reserved.
17 !
18 ! FVCOM has been developed by the joint UMASSD-WHOI research team. For
19 ! details of authorship and attribution of credit please see the FVCOM
20 ! technical manual or contact the MEDM group.
21 !
22 !
23 ! This file is part of FVCOM. For details, see http://fvcom.smast.umassd.edu
24 ! The full copyright notice is contained in the file COPYRIGHT located in the
25 ! root directory of the FVCOM code. This original header must be maintained
26 ! in all distributed versions.
27 !
28 ! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
29 ! AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO,
30 ! THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
31 ! PURPOSE ARE DISCLAIMED.
32 !
33 !/---------------------------------------------------------------------------/
34 ! CVS VERSION INFORMATION
35 ! $Id$
36 ! $Name$
37 ! $Revision$
38 !/===========================================================================/
39 
40 !==============================================================================!
41 ! LAGRANGIAN PARTICLE CLASS !
42 !==============================================================================!
43 
44 ! MODFIED BY DAVID FOR THE NEW INPUT: NEEDS TESTING
45 
46 ! ALL TIMES ARE NOW EITHER TYPE(TIME) OR IN DAYS!
47 
49  use mod_prec
50  use mod_time
51  use control, only : zerotime
52  USE mod_utils
53 
54  ! INTEL 9.1.X does not like initialized types
55  integer, parameter :: number_of_floats = 26
56  integer, parameter :: number_of_scalars = 5
57 
58 !!$ type particle
59 !!$ integer :: id !!global particle number
60 !!$ integer :: pid !!processor id number
61 !!$ integer :: elem !!element containing particle
62 !!$ integer :: group !!element group id
63 !!$ TYPE(TIME) :: tbeg !! START TIME
64 !!$ TYPE(TIME) :: tend !! END TIME
65 !!$ real(sp) :: x(3) !!particle position
66 !!$ real(sp) :: xn(3) !!last time step position
67 !!$ real(sp) :: el !!surface elevation at particle point
68 !!$ real(sp) :: h !!bathymetry at particle point
69 !!$ real(sp) :: u !!x-velocity at particle location
70 !!$ real(sp) :: v !!y-velocity at particle location
71 !!$ real(sp) :: w !!sigma-velocity at particle location
72 !!$ real(sp) :: zloc !!z position of particle [m]
73 !!$ real(sp) :: s(number_of_scalars) !!particle scalar
74 !!$ real(sp) :: chi(3,4) !!Runge-Kutta stage contributions
75 !!$ real(sp) :: deltat !!particle time step
76 !!$ real(sp) :: pathlength !!particle integrated pathlength [m]
77 !!$ logical :: found
78 !!$ end type particle
79 
80  ! USE THIS FOR INTEL 10.1+
81  type particle
82  integer :: id =0 !!global particle number
83  integer :: pid =0 !!processor id number
84  integer :: elem =0 !!element containing particle
85  integer :: group =0 !!element group id
86  TYPE(time) :: tbeg !! START TIME
87  TYPE(time) :: tend !! END TIME
88  real(sp) :: x(3) =0.0_sp !!particle position
89  real(sp) :: xn(3) =0.0_sp !!last time step position
90  real(sp) :: el =0._sp !!surface elevation at particle point
91  real(sp) :: h =0._sp !!bathymetry at particle point
92  real(sp) :: u =0._sp !!x-velocity at particle location
93  real(sp) :: v =0._sp !!y-velocity at particle location
94  real(sp) :: w =0._sp !!sigma-velocity at particle location
95  real(sp) :: zloc =0._sp !!z position of particle [m]
96  real(sp) :: s(number_of_scalars) =0._sp !!particle scalar
97  real(sp) :: chi(3,4) =0._sp !!Runge-Kutta stage contributions
98  real(sp) :: deltat =0._sp !!particle time step
99  real(sp) :: pathlength =0._sp !!particle integrated pathlength [m]
100  logical :: found = .false.
101  end type particle
102 
103  INTEGER mpi_particle
104 
105  interface operator (<) !for sorting, insert
106  module procedure less_than_particle ; end interface
107  interface operator (==) !for sorting or delete
108  module procedure equal_to_particle ; end interface
109 
110  contains
111 
112  function less_than_particle (particle_1,particle_2) result(Boolean)
113  type(particle), intent(in) :: particle_1
114  type(particle), intent(in) :: particle_2
115  logical :: boolean
116  boolean = particle_1%id < particle_2%id
117  end function less_than_particle
118 
119  function equal_to_particle (particle_1,particle_2) result(Boolean)
120  type(particle), intent(in) :: particle_1
121  type(particle), intent(in) :: particle_2
122  logical :: boolean
123  boolean = particle_1%id == particle_2%id
124  end function equal_to_particle
125 
126  subroutine screen_write(p,hprint)
127  type(particle), intent(in) :: p
128  logical, intent(in) :: hprint
129  if(hprint)then
130  write(ipt,*)'id group x y z elem pid'
131  endif
132  write(ipt,'(2I5,3F10.2,I8,2I5)')p%id,p%group,p%x,p%elem,p%pid
133  end subroutine screen_write
134 
135  subroutine shift_pos(p)
136  type(particle), intent(inout) :: p
137  p%xn = p%x
138  end subroutine shift_pos
139 
140 !!$ subroutine set_pathlength(p)
141 !!$ type(particle), intent(inout) :: p
142 !!$ p%pathlength = sqrt( (p%x-p%xn)**2)
143 !!$ end subroutine set_pathlength
144 
145  !dump particle info to screen
146  subroutine particle_print(p)
147  type(particle), intent(in) :: p
148 ! type(particle), pointer :: p
149  write(ipt,*)
150  write(ipt,*)'id: ',p%id
151  write(ipt,*)'processor: ',p%pid
152  write(ipt,*)'element: ',p%elem
153  write(ipt,*)'group: ',p%group
154  write(ipt,*)'tbeg: ',p%tbeg
155  write(ipt,*)'tend: ',p%tend
156  write(ipt,*)'x: ',p%x
157  write(ipt,*)'xn: ',p%xn
158  write(ipt,*)'el: ',p%el
159  write(ipt,*)'h: ',p%h
160  write(ipt,*)'u: ',p%u
161  write(ipt,*)'v: ',p%v
162  write(ipt,*)'w: ',p%w
163  write(ipt,*)'zloc: ',p%zloc
164  write(ipt,*)'s: ',p%s
165  write(ipt,*)'chi1 ',p%chi(:,1)
166  write(ipt,*)'chi2 ',p%chi(:,2)
167  write(ipt,*)'chi3 ',p%chi(:,3)
168  write(ipt,*)'chi4 ',p%chi(:,4)
169  write(ipt,*)'deltat ',p%deltat
170  write(ipt,*)'pathlength ',p%pathlength
171  end subroutine particle_print
172 
173  !initialize values
174  subroutine zero_out(p)
175  implicit none
176  type(particle), intent(inout) :: p
177 
178  ! SET VALUES
179  p%id = 0
180  p%pid = 0
181  p%group = 0
182  p%elem = 0
183  p%tbeg = zerotime
184  p%tend = zerotime
185  p%x = 0.0_sp
186  p%xn = 0.0_sp
187  p%el = 0.0_sp
188  p%u = 0.0_sp
189  p%v = 0.0_sp
190  p%w = 0.0_sp
191  p%zloc = 0.0_sp
192  p%s = 0.0_sp
193  p%chi = 0.0_sp
194  p%deltat = 0.0
195  p%pathlength = 0.0
196  p%found = .false.
197  end subroutine zero_out
198 
199  SUBROUTINE new_particles(P,i)
200  implicit none
201  integer, intent(in) :: i
202  type(particle), allocatable, intent(inout) :: p(:)
203  integer :: status, J
204 
205  allocate(p(i),stat=status)
206  if(status/=0) CALL fatal_error("NEW_PARTICLE: COULD NOT ALLOCATE?")
207 
208 
209  DO j=1,i
210  call zero_out(p(j))
211  END DO
212 
213  END SUBROUTINE new_particles
214 
215 
216  SUBROUTINE create_mpi_particle
218  IMPLICIT NONE
219  END SUBROUTINE create_mpi_particle
220 
221 end module particle_class
subroutine new_particles(P, i)
Definition: particle.f90:200
subroutine create_mpi_particle
Definition: particle.f90:217
subroutine screen_write(p, hprint)
Definition: particle.f90:127
type(time) zerotime
Definition: mod_main.f90:830
logical function less_than_particle(particle_1, particle_2)
Definition: particle.f90:113
subroutine shift_pos(p)
Definition: particle.f90:136
integer, parameter number_of_floats
Definition: particle.f90:55
integer mpi_particle
Definition: particle.f90:103
integer, parameter sp
Definition: mod_prec.f90:48
integer, parameter number_of_scalars
Definition: particle.f90:56
subroutine fatal_error(ER1, ER2, ER3, ER4)
Definition: mod_utils.f90:230
subroutine zero_out(p)
Definition: particle.f90:175
logical function equal_to_particle(particle_1, particle_2)
Definition: particle.f90:120
subroutine particle_print(p)
Definition: particle.f90:147