~ubuntu-branches/ubuntu/precise/code-saturne/precise

« back to all changes in this revision

Viewing changes to users/base/ustsns.f90

  • Committer: Package Import Robot
  • Author(s): Sylvestre Ledru
  • Date: 2011-11-24 00:00:08 UTC
  • mfrom: (6.1.9 sid)
  • Revision ID: package-import@ubuntu.com-20111124000008-2vo99e38267942q5
Tags: 2.1.0-3
Install a missing file

Show diffs side-by-side

added added

removed removed

Lines of Context:
2
2
 
3
3
!VERS
4
4
 
5
 
 
6
 
!     This file is part of the Code_Saturne Kernel, element of the
7
 
!     Code_Saturne CFD tool.
8
 
 
9
 
!     Copyright (C) 1998-2009 EDF S.A., France
10
 
 
11
 
!     contact: saturne-support@edf.fr
12
 
 
13
 
!     The Code_Saturne Kernel is free software; you can redistribute it
14
 
!     and/or modify it under the terms of the GNU General Public License
15
 
!     as published by the Free Software Foundation; either version 2 of
16
 
!     the License, or (at your option) any later version.
17
 
 
18
 
!     The Code_Saturne Kernel is distributed in the hope that it will be
19
 
!     useful, but WITHOUT ANY WARRANTY; without even the implied warranty
20
 
!     of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
21
 
!     GNU General Public License for more details.
22
 
 
23
 
!     You should have received a copy of the GNU General Public License
24
 
!     along with the Code_Saturne Kernel; if not, write to the
25
 
!     Free Software Foundation, Inc.,
26
 
!     51 Franklin St, Fifth Floor,
27
 
!     Boston, MA  02110-1301  USA
 
5
! This file is part of Code_Saturne, a general-purpose CFD tool.
 
6
!
 
7
! Copyright (C) 1998-2011 EDF S.A.
 
8
!
 
9
! This program is free software; you can redistribute it and/or modify it under
 
10
! the terms of the GNU General Public License as published by the Free Software
 
11
! Foundation; either version 2 of the License, or (at your option) any later
 
12
! version.
 
13
!
 
14
! This program is distributed in the hope that it will be useful, but WITHOUT
 
15
! ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
 
16
! FOR A PARTICULAR PURPOSE.  See the GNU General Public License for more
 
17
! details.
 
18
!
 
19
! You should have received a copy of the GNU General Public License along with
 
20
! this program; if not, write to the Free Software Foundation, Inc., 51 Franklin
 
21
! Street, Fifth Floor, Boston, MA 02110-1301, USA.
28
22
 
29
23
!-------------------------------------------------------------------------------
30
24
 
31
25
subroutine ustsns &
32
26
!================
33
27
 
34
 
 ( idbia0 , idbra0 ,                                              &
35
 
   ndim   , ncelet , ncel   , nfac   , nfabor , nfml   , nprfml , &
36
 
   nnod   , lndfac , lndfbr , ncelbr ,                            &
37
 
   nvar   , nscal  , nphas  , ncepdp , ncesmp ,                   &
38
 
   nideve , nrdeve , nituse , nrtuse ,                            &
39
 
   ivar   , iphas  ,                                              &
40
 
   ifacel , ifabor , ifmfbr , ifmcel , iprfml , maxelt , lstelt , &
41
 
   ipnfac , nodfac , ipnfbr , nodfbr ,                            &
 
28
 ( nvar   , nscal  , ncepdp , ncesmp ,                            &
 
29
   ivar   ,                                                       &
42
30
   icepdc , icetsm , itypsm ,                                     &
43
 
   idevel , ituser , ia     ,                                     &
44
 
   xyzcen , surfac , surfbo , cdgfac , cdgfbo , xyznod , volume , &
45
31
   dt     , rtpa   , propce , propfa , propfb ,                   &
46
32
   coefa  , coefb  , ckupdc , smacel ,                            &
47
 
   crvexp , crvimp ,                                              &
48
 
   dam    , xam    ,                                              &
49
 
   w1     , w2     , w3     , w4     , w5     , w6     ,          &
50
 
   rdevel , rtuser , ra     )
 
33
   crvexp , crvimp )
51
34
 
52
35
!===============================================================================
53
36
! Purpose:
63
46
! -----
64
47
! The routine is called for each velocity component. It is therefore necessary
65
48
! to test the value of the variable ivar to separate the treatments of the
66
 
! components iu(iphas), iv(iphas) or iw(iphas).
 
49
! components iu, iv or iw.
67
50
!
68
51
! The additional source term is decomposed into an explicit part (crvexp) and
69
52
! an implicit part (crvimp) that must be provided here.
101
84
!__________________.____._____.________________________________________________.
102
85
! name             !type!mode ! role                                           !
103
86
!__________________!____!_____!________________________________________________!
104
 
! idbia0           ! i  ! <-- ! number of first free position in ia            !
105
 
! idbra0           ! i  ! <-- ! number of first free position in ra            !
106
 
! ndim             ! i  ! <-- ! spatial dimension                              !
107
 
! ncelet           ! i  ! <-- ! number of extended (real + ghost) cells        !
108
 
! ncel             ! i  ! <-- ! number of cells                                !
109
 
! nfac             ! i  ! <-- ! number of interior faces                       !
110
 
! nfabor           ! i  ! <-- ! number of boundary faces                       !
111
 
! nfml             ! i  ! <-- ! number of families (group classes)             !
112
 
! nprfml           ! i  ! <-- ! number of properties per family (group class)  !
113
 
! nnod             ! i  ! <-- ! number of vertices                             !
114
 
! lndfac           ! i  ! <-- ! size of nodfac indexed array                   !
115
 
! lndfbr           ! i  ! <-- ! size of nodfbr indexed array                   !
116
 
! ncelbr           ! i  ! <-- ! number of cells with faces on boundary         !
117
87
! nvar             ! i  ! <-- ! total number of variables                      !
118
88
! nscal            ! i  ! <-- ! total number of scalars                        !
119
 
! nphas            ! i  ! <-- ! number of phases                               !
120
89
! ncepdp           ! i  ! <-- ! number of cells with head loss terms           !
121
90
! ncssmp           ! i  ! <-- ! number of cells with mass source terms         !
122
 
! nideve, nrdeve   ! i  ! <-- ! sizes of idevel and rdevel arrays              !
123
 
! nituse, nrtuse   ! i  ! <-- ! sizes of ituser and rtuser arrays              !
124
91
! ivar             ! i  ! <-- ! index number of the current variable           !
125
 
! iphas            ! i  ! <-- ! index number of the current phase              !
126
 
! ifacel(2, nfac)  ! ia ! <-- ! interior faces -> cells connectivity           !
127
 
! ifabor(nfabor)   ! ia ! <-- ! boundary faces -> cells connectivity           !
128
 
! ifmfbr(nfabor)   ! ia ! <-- ! boundary face family numbers                   !
129
 
! ifmcel(ncelet)   ! ia ! <-- ! cell family numbers                            !
130
 
! iprfml           ! ia ! <-- ! property numbers per family                    !
131
 
!  (nfml, nprfml)  !    !     !                                                !
132
 
! maxelt           ! i  ! <-- ! max number of cells and faces (int/boundary)   !
133
 
! lstelt(maxelt)   ! ia ! --- ! work array                                     !
134
 
! ipnfac(nfac+1)   ! ia ! <-- ! interior faces -> vertices index (optional)    !
135
 
! nodfac(lndfac)   ! ia ! <-- ! interior faces -> vertices list (optional)     !
136
 
! ipnfbr(nfabor+1) ! ia ! <-- ! boundary faces -> vertices index (optional)    !
137
 
! nodfbr(lndfbr)   ! ia ! <-- ! boundary faces -> vertices list (optional)     !
138
92
! icepdc(ncepdp)   ! ia ! <-- ! index number of cells with head loss terms     !
139
93
! icetsm(ncesmp)   ! ia ! <-- ! index number of cells with mass source terms   !
140
94
! itypsm           ! ia ! <-- ! type of mass source term for each variable     !
141
95
!  (ncesmp,nvar)   !    !     !  (see ustsma)                                  !
142
 
! idevel(nideve)   ! ia ! <-- ! integer work array for temporary developpement !
143
 
! ituser(nituse    ! ia ! <-- ! user-reserved integer work array               !
144
 
! ia(*)            ! ia ! --- ! main integer work array                        !
145
 
! xyzcen           ! ra ! <-- ! cell centers                                   !
146
 
!  (ndim, ncelet)  !    !     !                                                !
147
 
! surfac           ! ra ! <-- ! interior faces surface vectors                 !
148
 
!  (ndim, nfac)    !    !     !                                                !
149
 
! surfbo           ! ra ! <-- ! boundary faces surface vectors                 !
150
 
!  (ndim, nfavor)  !    !     !                                                !
151
 
! cdgfac           ! ra ! <-- ! interior faces centers of gravity              !
152
 
!  (ndim, nfac)    !    !     !                                                !
153
 
! cdgfbo           ! ra ! <-- ! boundary faces centers of gravity              !
154
 
!  (ndim, nfabor)  !    !     !                                                !
155
 
! xyznod           ! ra ! <-- ! vertex coordinates (optional)                  !
156
 
!  (ndim, nnod)    !    !     !                                                !
157
 
! volume(ncelet)   ! ra ! <-- ! cell volumes                                   !
158
96
! dt(ncelet)       ! ra ! <-- ! time step (per cell)                           !
159
97
! rtpa             ! ra ! <-- ! calculated variables at cell centers           !
160
98
!  (ncelet, *)     !    !     !  (preceding time steps)                        !
168
106
!  (ncesmp,nvar)   !    !     !  source terms or mass rate (see ustsma)        !
169
107
! crvexp           ! ra ! --> ! explicit part of the source term               !
170
108
! crvimp           ! ra ! --> ! implicit part of the source term               !
171
 
! dam(ncelet)      ! ra ! --- ! work array                                     !
172
 
! xam(nfac,2)      ! ra ! --- ! work array                                     !
173
 
! w1,2,3,4,5,6     ! ra ! --- ! work arrays                                    !
174
 
!  (ncelet)        !    !     !  (computation of pressure gradient)            !
175
 
! rdevel(nrdeve)   ! ra ! <-> ! real work array for temporary developpement    !
176
 
! rtuser(nituse    ! ra ! <-- ! user-reserved real work array                  !
177
 
! ra(*)            ! ra ! --- ! main real work array                           !
178
109
!__________________!____!_____!________________________________________________!
179
110
 
180
111
!     Type: i (integer), r (real), s (string), a (array), l (logical),
182
113
!     mode: <-- input, --> output, <-> modifies data, --- work array
183
114
!===============================================================================
184
115
 
 
116
!===============================================================================
 
117
! Module files
 
118
!===============================================================================
 
119
 
 
120
use paramx
 
121
use numvar
 
122
use entsor
 
123
use optcal
 
124
use cstphy
 
125
use parall
 
126
use period
 
127
use mesh
 
128
 
 
129
!===============================================================================
 
130
 
185
131
implicit none
186
132
 
187
 
!===============================================================================
188
 
! Common blocks
189
 
!===============================================================================
190
 
 
191
 
include "paramx.h"
192
 
include "pointe.h"
193
 
include "numvar.h"
194
 
include "entsor.h"
195
 
include "optcal.h"
196
 
include "cstphy.h"
197
 
include "parall.h"
198
 
include "period.h"
199
 
 
200
 
!===============================================================================
201
 
 
202
133
! Arguments
203
134
 
204
 
integer          idbia0 , idbra0
205
 
integer          ndim   , ncelet , ncel   , nfac   , nfabor
206
 
integer          nfml   , nprfml
207
 
integer          nnod   , lndfac , lndfbr , ncelbr
208
 
integer          nvar   , nscal  , nphas
 
135
integer          nvar   , nscal
209
136
integer          ncepdp , ncesmp
210
 
integer          nideve , nrdeve , nituse , nrtuse
211
 
integer          ivar   , iphas
 
137
integer          ivar
212
138
 
213
 
integer          ifacel(2,nfac) , ifabor(nfabor)
214
 
integer          ifmfbr(nfabor) , ifmcel(ncelet)
215
 
integer          iprfml(nfml,nprfml), maxelt, lstelt(maxelt)
216
 
integer          ipnfac(nfac+1), nodfac(lndfac)
217
 
integer          ipnfbr(nfabor+1), nodfbr(lndfbr)
218
139
integer          icepdc(ncepdp)
219
140
integer          icetsm(ncesmp), itypsm(ncesmp,nvar)
220
 
integer          idevel(nideve), ituser(nituse), ia(*)
221
141
 
222
 
double precision xyzcen(ndim,ncelet)
223
 
double precision surfac(ndim,nfac), surfbo(ndim,nfabor)
224
 
double precision cdgfac(ndim,nfac), cdgfbo(ndim,nfabor)
225
 
double precision xyznod(ndim,nnod), volume(ncelet)
226
142
double precision dt(ncelet), rtpa(ncelet,*)
227
143
double precision propce(ncelet,*)
228
144
double precision propfa(nfac,*), propfb(nfabor,*)
229
145
double precision coefa(nfabor,*), coefb(nfabor,*)
230
146
double precision ckupdc(ncepdp,6), smacel(ncesmp,nvar)
231
147
double precision crvexp(ncelet), crvimp(ncelet)
232
 
double precision dam(ncelet ),xam(nfac ,2)
233
 
double precision w1(ncelet),w2(ncelet),w3(ncelet)
234
 
double precision w4(ncelet),w5(ncelet),w6(ncelet)
235
 
double precision rdevel(nrdeve), rtuser(nrtuse), ra(*)
236
148
 
237
149
! Local variables
238
150
 
239
151
character*80     chaine
240
 
integer          idebia, idebra
241
152
integer          iel, ipcrom, ipp, iutile
242
153
double precision ckp, qdm
243
154
 
 
155
integer, allocatable, dimension(:) :: lstelt
 
156
 
244
157
!===============================================================================
245
158
 
246
159
! TEST_TO_REMOVE_FOR_USE_OF_SUBROUTINE_START
256
169
! 1. Initialization
257
170
!===============================================================================
258
171
 
259
 
idebia = idbia0
260
 
idebra = idbra0
 
172
! Allocate a temporary array for cells selection
 
173
allocate(lstelt(ncel))
 
174
 
261
175
 
262
176
ipp    = ipprtp(ivar)
263
177
 
266
180
  write(nfecra,1000) chaine(1:8)
267
181
endif
268
182
 
269
 
ipcrom = ipproc(irom  (iphas))
 
183
ipcrom = ipproc(irom  )
270
184
 
271
185
 
272
186
!===============================================================================
303
217
 
304
218
! ----------------------------------------------
305
219
 
306
 
iphas=1
307
 
if (ivar.eq.iu(1)) then
 
220
if (ivar.eq.iu) then
308
221
 
309
222
  ckp  = 10.d0
310
223
  qdm  = 100.d0
330
243
! End
331
244
!----
332
245
 
 
246
! Deallocate the temporary array
 
247
deallocate(lstelt)
 
248
 
333
249
return
334
 
 
335
250
end subroutine