File: cscpce.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 (149 lines) | stat: -rw-r--r-- 4,240 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
!-------------------------------------------------------------------------------

! 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.

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

!> \file cscpce.f90
!> \brief Preparation of sending velocity variables for coupling between
!> two instances of Code_Saturne via boundary faces.
!> Received indformation will be transformed into boundary condition
!> in subroutine \ref csc2cl.
!------------------------------------------------------------------------------

!------------------------------------------------------------------------------
! Arguments
!------------------------------------------------------------------------------
!   mode          name          role
!------------------------------------------------------------------------------
!> \param[in]     nptdis
!> \param[in]     ivar          variable number
!> \param[in]     locpts
!> \param[in]     vela          variable value at time step beginning
!> \param[in]     coefav
!> \param[in]     coefbv
!> \param[in]     coopts
!> \param[out]    rvdis
!______________________________________________________________________________

subroutine cscpce &
 ( nptdis , ivar   ,                                              &
   locpts ,                                                       &
   vela   ,                                                       &
   coefav , coefbv ,                                              &
   coopts , rvdis  )

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

use paramx
use pointe
use numvar
use optcal
use cstphy
use cstnum
use entsor
use parall
use period
use cplsat
use mesh

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

implicit none

! Arguments

integer          ivar
integer          nptdis

integer          locpts(nptdis)

double precision coopts(3,nptdis), rvdis(3,nptdis)
double precision coefav(3  ,nfabor)
double precision coefbv(3,3,nfabor)
double precision vela(3,ncelet)

! Local variables

integer          ipt    , iel    , isou   , f_id
integer          inc    , iccocg , nswrgp
integer          iwarnp , imligp

double precision epsrgp , climgp
double precision dx     , dy     , dz

double precision, dimension(:,:,:), allocatable :: gradv

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

! Allocate a temporary array
allocate(gradv(3,3,ncelet))

inc    = 1
iccocg = 1
nswrgp = nswrgr(ivar)
imligp = imligr(ivar)
iwarnp = iwarni(ivar)
epsrgp = epsrgr(ivar)
climgp = climgr(ivar)

if (ivar.le.0) then
  f_id = -1
else
  f_id = ivarfl(ivar)
endif

call cgdvec &
( f_id   , imrgra , inc    , nswrgp , iwarnp , imligp ,          &
  epsrgp , climgp ,                                              &
  coefav , coefbv , vela   ,                                     &
  gradv)

! --- Interpolation

do ipt = 1, nptdis

  iel = locpts(ipt)

  dx = coopts(1,ipt) - xyzcen(1,iel)
  dy = coopts(2,ipt) - xyzcen(2,iel)
  dz = coopts(3,ipt) - xyzcen(3,iel)

  do isou = 1, 3
    rvdis(isou,ipt) = vela(isou,iel) + gradv(1,isou,iel)*dx       &
                                     + gradv(2,isou,iel)*dy       &
                                     + gradv(3,isou,iel)*dz
  enddo

enddo

! Free memory
deallocate(gradv)

!--------
! FORMATS
!--------
!----
! FIN
!----

return
end subroutine