File: strdep.f90

package info (click to toggle)
code-saturne 4.3.3%2Brepack-1
  • links: PTS, VCS
  • area: main
  • in suites: stretch
  • size: 77,992 kB
  • sloc: ansic: 281,257; f90: 122,305; python: 56,490; makefile: 3,915; xml: 3,285; cpp: 3,183; sh: 1,139; lex: 176; yacc: 101; sed: 16
file content (391 lines) | stat: -rw-r--r-- 12,969 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
!-------------------------------------------------------------------------------

! This file is part of Code_Saturne, a general-purpose CFD tool.
!
! Copyright (C) 1998-2016 EDF S.A.
!
! This program is free software; you can redistribute it and/or modify it under
! the terms of the GNU General Public License as published by the Free Software
! Foundation; either version 2 of the License, or (at your option) any later
! version.
!
! This program is distributed in the hope that it will be useful, but WITHOUT
! ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
! FOR A PARTICULAR PURPOSE.  See the GNU General Public License for more
! details.
!
! You should have received a copy of the GNU General Public License along with
! this program; if not, write to the Free Software Foundation, Inc., 51 Franklin
! Street, Fifth Floor, Boston, MA 02110-1301, USA.

!-------------------------------------------------------------------------------

subroutine strdep &
!================

 ( itrale , italim , itrfin ,                                     &
   nvar   ,                                                       &
   dt     ,                                                       &
   flmalf , flmalb , cofale , xprale )

!===============================================================================
! FONCTION :
! ----------

! DEPLACEMENT DES STRUCTURES MOBILES EN ALE EN COUPLAGE INTERNE

!-------------------------------------------------------------------------------
! Arguments
!__________________.____._____.________________________________________________.
! name             !type!mode ! role                                           !
!__________________!____!_____!________________________________________________!
! itrale           ! e  ! <-- ! numero d'iteration pour l'ale                  !
! italim           ! e  ! <-- ! numero d'iteration couplage implicite          !
! itrfin           ! e  ! <-- ! indicateur de derniere iteration de            !
!                  !    !     !                    couplage implicite          !
! nvar             ! i  ! <-- ! total number of variables                      !
! dt(ncelet)       ! ra ! <-- ! time step (per cell)                           !
! flmalf(nfac)     ! tr ! --> ! sauvegarde du flux de masse faces int          !
! flmalb(nfabor    ! tr ! --> ! sauvegarde du flux de masse faces brd          !
! cofale           ! tr ! --> ! sauvegarde des cl de p et u                    !
!    (nfabor,8)    !    !     !                                                !
! xprale(ncelet    ! tr ! --> ! sauvegarde de la pression, si nterup           !
!                  !    !     !    est >1                                      !
!__________________!____!_____!________________________________________________!

!     TYPE : E (ENTIER), R (REEL), A (ALPHANUMERIQUE), T (TABLEAU)
!            L (LOGIQUE)   .. ET TYPES COMPOSES (EX : TR TABLEAU REEL)
!     MODE : <-- donnee, --> resultat, <-> Donnee modifiee
!            --- tableau de travail
!===============================================================================

!===============================================================================
! Module files
!===============================================================================

use paramx
use ihmpre
use cstphy
use numvar
use optcal
use entsor
use pointe
use albase
use alstru
use alaste
use parall
use period
use mesh
use field

!===============================================================================

implicit none

! Arguments

integer          itrale , italim , itrfin
integer          nvar

double precision dt(ncelet)
double precision flmalf(nfac), flmalb(nfabor), xprale(ncelet)
double precision cofale(nfabor,11)

! Local variables

integer          istr, ii, iel, ifac, ntab
integer          iflmas, iflmab
integer          indast
integer          icvext, icvint, icved
integer          f_dim

double precision delta

double precision, allocatable, dimension(:,:) :: forast
double precision, dimension(:), pointer :: imasfl, bmasfl
double precision, dimension(:,:), pointer :: forbr
double precision, dimension(:,:), pointer :: coefau
double precision, dimension(:,:,:), pointer :: coefbu
double precision, dimension(:), pointer :: coefap, coefbp
double precision, dimension(:), pointer :: cvar_var, cvara_var
double precision, dimension(:,:), pointer :: cvar_varv, cvara_varv

!===============================================================================

!===============================================================================
! 1. Initialization
!===============================================================================

call field_get_key_int(ivarfl(iu), kimasf, iflmas)
call field_get_key_int(ivarfl(iu), kbmasf, iflmab)
call field_get_val_s(iflmas, imasfl)
call field_get_val_s(iflmab, bmasfl)

call field_get_val_v(iforbr, forbr)

call field_get_coefa_v(ivarfl(iu), coefau)
call field_get_coefb_v(ivarfl(iu), coefbu)

call field_get_coefa_s(ivarfl(ipr), coefap)
call field_get_coefb_s(ivarfl(ipr), coefbp)

!===============================================================================
! 2. Computue forces on the structures
!===============================================================================

do istr = 1, nbstru
  do ii = 1, ndim
    forsta(ii,istr) = forstr(ii,istr)
!-a tester          forsta(ii,istr) = forstP(ii,istr)
    forstr(ii,istr) = 0.d0
  enddo
enddo

! Allocate a temporary array
allocate(forast(3,nbfast))

indast = 0
do ifac = 1, nfabor
  istr = idfstr(ifac)
  if (istr.gt.0) then
    do ii = 1, 3
      forstr(ii,istr) = forstr(ii,istr) + forbr(ii,ifac)
    enddo
  else if (istr.lt.0) then
    indast = indast + 1
    do ii = 1, 3
      forast(ii,indast) = asddlf(ii,-istr)*forbr(ii,ifac)
    enddo
  endif
enddo

if (irangp.ge.0) then
  ntab = ndim*nbstru
  call parrsm(ntab,forstr)
endif

!     Calcul de l'effort envoye au structures internes
do istr = 1, nbstru
  do ii = 1, ndim
    forstp(ii,istr) = cfopre*forstr(ii,istr)+                     &
         (1.d0-cfopre)*forsta(ii,istr)
  enddo
enddo

!     Envoi de l'effort applique aux structures externes
if (nbaste.gt.0) then
  call astfor(ntcast, nbfast, forast)
  !==========
endif

! Free memory
deallocate(forast)

!===============================================================================
! 3. Structure characteristic given by the user
!===============================================================================

if (nbstru.gt.0) then

  ! - Interface Code_Saturne
  !   ======================

  if (iihmpr.eq.1) then

    call uistr2 &
    !==========
 ( xmstru, xcstru, xkstru,     &
   forstp,                     &
   dtref , ttcabs, ntcabs   )

  endif

  call usstr2                                                     &
  !==========
 ( nbstru ,                                                       &
   idfstr ,                                                       &
   dt     ,                                                       &
   xmstru , xcstru , xkstru , xstreq , xstr   , xpstr  , forstp , &
   dtstr  )

endif

! If the fluid is initailizing, we do read structures
if (itrale.le.nalinf) then
  itrfin = -1
  return
endif

!===============================================================================
! 4. Internal structure displacement
!===============================================================================

do istr = 1, nbstru

  call newmrk                                                     &
  !==========
 ( istr  , alpnmk  , betnmk          , gamnmk          ,          &
   xmstru(1,1,istr), xcstru(1,1,istr), xkstru(1,1,istr),          &
   xstreq(1,istr)  ,                                              &
   xstr(1,istr)    , xpstr(1,istr)   , xppstr(1,istr)  ,          &
   xsta(1,istr)    , xpsta(1,istr)   , xppsta(1,istr)  ,          &
   forstp(1,istr)  , forsta(1,istr)  , dtstr(istr)     )

enddo

!===============================================================================
! 5. Convergence test
!===============================================================================

icvext = 0
icvint = 0
icved  = 0

delta = 0.d0
do istr = 1, nbstru
  do ii = 1, 3
    delta = delta + (xstr(ii,istr)-xstp(ii,istr))**2
  enddo
enddo
if (nbstru.gt.0) then
  delta = sqrt(delta)/almax/nbstru
  if (delta.lt.epalim) icvint = 1
endif

if (nbaste.gt.0) call astcv1(ntcast, icvext)

if (nbstru.gt.0.and.nbaste.gt.0) then
   icved = icvext*icvint
elseif (nbstru.gt.0.and.nbaste.eq.0) then
   icved = icvint
elseif (nbaste.gt.0.and.nbstru.eq.0) then
   icved = icvext
endif

if (iwarni(iuma).ge.2) write(nfecra,1000) italim, delta

!     si convergence
if (icved.eq.1) then
  if (itrfin.eq.1) then
!       si ITRFIN=1 on sort
    if (iwarni(iuma).ge.1) write(nfecra,1001) italim, delta
    itrfin = -1
  else
!       sinon on refait une derniere iteration pour SYRTHES/T1D/rayonnement
!        et on remet ICVED a 0 pour que Code_Aster refasse une iteration aussi
    itrfin = 1
    icved = 0
  endif
elseif (itrfin.eq.0 .and. italim.eq.nalimx-1) then
!       ce sera la derniere iteration
  itrfin = 1
elseif (italim.eq.nalimx) then
!       on a forcement ITRFIN=1 et on sort
  if (nalimx.gt.1) write(nfecra,1100) italim, delta
  itrfin = -1
!       On met ICVED a 1 pour que Code_Aster s'arrete lui aussi
  icved = 1
endif

!     On renvoie l'indicateur de convergence final a Code_Aster
call astcv2(ntcast, icved)
!==========

!===============================================================================
! 6. Re set previous values if required
!===============================================================================

! If NTERUP .GT. 1, values at previous time step have been modified after NAVSTO
! We must then go back to a previous value
if (itrfin.ne.-1) then
  do ii = 1, nvar
    call field_get_dim (ivarfl(ii), f_dim)

    ! Fields of dimension one
    if (f_dim.eq.1) then
      call field_get_val_s(ivarfl(ii), cvar_var)
      call field_get_val_prev_s(ivarfl(ii), cvara_var)
      if (ii.eq.ipr .and. nterup.gt.1) then
        do iel = 1, ncelet
          cvara_var(iel) = xprale(iel)
        enddo
      endif
      do iel = 1, ncelet
        cvar_var(iel) = cvara_var(iel)
      enddo

    ! Vector fields
    else if (f_dim.eq.3) then
      call field_get_val_v(ivarfl(ii), cvar_varv)
      call field_get_val_prev_v(ivarfl(ii), cvara_varv)
      do iel = 1, ncelet
        cvar_varv(1, iel) = cvara_varv(1, iel)
        cvar_varv(2, iel) = cvara_varv(2, iel)
        cvar_varv(3, iel) = cvara_varv(3, iel)
      enddo
    else
      call csexit(1)
    endif
  enddo
  do ifac = 1, nfac
    imasfl(ifac) = flmalf(ifac)
  enddo
  do ifac = 1, nfabor
    bmasfl(ifac) = flmalb(ifac)
    coefap(ifac) = cofale(ifac,1)
    coefau(1, ifac) = cofale(ifac,2)
    coefau(2, ifac) = cofale(ifac,3)
    coefau(3, ifac) = cofale(ifac,4)
    coefbp(ifac) = cofale(ifac,5)
    coefbu(1, 1, ifac) = cofale(ifac,6)
    coefbu(2, 2, ifac) = cofale(ifac,7)
    coefbu(3, 3, ifac) = cofale(ifac,8)
    coefbu(1, 2, ifac) = cofale(ifac,9)
    coefbu(2, 3, ifac) = cofale(ifac,10)
    coefbu(1, 3, ifac) = cofale(ifac,11)
    ! the coefficient B is supposed to be symmetric
    coefbu(2, 1, ifac) = cofale(ifac,9)
    coefbu(3, 2, ifac) = cofale(ifac,10)
    coefbu(3, 1, ifac) = cofale(ifac,11)
  enddo
endif

!--------
! Formats
!--------

#if defined(_CS_LANG_FR)

 1000 format (                                                          &
 '            ALE IMPLICITE : ITER=',I5,' DERIVE=',E12.5     )
 1001 format (                                                          &
 'CONVERGENCE ALE IMPLICITE : ITER=',I5,' DERIVE=',E12.5     )
 1100 format (                                                          &
'@                                                            ',/,&
'@ @@ ATTENTION : COUPLAGE IMPLICITE ALE                      ',/,&
'@    =========                                               ',/,&
'@  Nombre d''iterations maximal ',I10   ,' atteint           ',/,&
'@  Derive normee :',E12.5                                     ,/,&
'@                                                            '  )

#else

 1000 format (                                                          &
 '            IMPLICIT ALE: ITER=',I5,' DERIVE=',E12.5     )
 1001 format (                                                          &
 'CONVERGENCE IMPLICIT ALE: ITER=',I5,' DERIVE=',E12.5     )
 1100 format (                                                          &
'@                                                            ',/,&
'@ @@ WARNING: IMPLICIT ALE                                   ',/,&
'@    ========                                                ',/,&
'@  Maximum number of iterations ',I10   ,' reached           ',/,&
'@  Normed derive :',E12.5                                     ,/,&
'@                                                            '  )

#endif

!----
! End
!----

end subroutine