xref: /phasta/phSolver/common/qpbc.f (revision 595995161822a203c8467e0e4a253d7bd7d6df32)
1*59599516SKenneth E. Jansen      subroutine qpbc( rmass, qres, iBC, iper, ilwork )
2*59599516SKenneth E. Jansenc---------------------------------------------------------------------
3*59599516SKenneth E. Jansenc
4*59599516SKenneth E. Jansenc This routine satisfies the periodic boundary conditions
5*59599516SKenneth E. Jansenc on the diffusive flux residual and mass matrix
6*59599516SKenneth E. Jansenc
7*59599516SKenneth E. Jansenc input:
8*59599516SKenneth E. Jansenc     rmass   (nshg)              : mass matrix
9*59599516SKenneth E. Jansenc     qres    (nshg,(nflow-1)*nsd) : diffusive flux vector
10*59599516SKenneth E. Jansenc
11*59599516SKenneth E. Jansenc output: modified qres and rmass
12*59599516SKenneth E. Jansenc---------------------------------------------------------------------
13*59599516SKenneth E. Jansen      include "common.h"
14*59599516SKenneth E. Jansen
15*59599516SKenneth E. Jansen      dimension rmass(nshg), qres(nshg,idflx),
16*59599516SKenneth E. Jansen     &          iBC(nshg), iper(nshg),uv(nshg,2),
17*59599516SKenneth E. Jansen     &          tmpvec(nshg,4), tmp(nshg)
18*59599516SKenneth E. Jansenc
19*59599516SKenneth E. Jansen      if(iabc==1) then   !are there any axisym bc's
20*59599516SKenneth E. Jansen      do i=1,idflx/nsd
21*59599516SKenneth E. Jansen         do j=1,2
22*59599516SKenneth E. Jansen            istrt=j+(i-1)*(nflow-1)
23*59599516SKenneth E. Jansen            uv(:,j)=qres(:,istrt)
24*59599516SKenneth E. Jansen         enddo
25*59599516SKenneth E. Jansen         call rotabc(uv, iBC, 'in ')
26*59599516SKenneth E. Jansen         do j=1,2
27*59599516SKenneth E. Jansen            istrt=j+(i-1)*(nflow-1)
28*59599516SKenneth E. Jansen            qres(:,istrt)=uv(:,j)
29*59599516SKenneth E. Jansen         enddo
30*59599516SKenneth E. Jansen      enddo
31*59599516SKenneth E. Jansen      endif
32*59599516SKenneth E. Jansenc
33*59599516SKenneth E. Jansenc
34*59599516SKenneth E. Jansenc.... compute qi for node A, i.e., qres <-- qres/rmass
35*59599516SKenneth E. Jansenc
36*59599516SKenneth E. Jansen       if (numpe > 1) then
37*59599516SKenneth E. Jansen          call commu (qres  , ilwork, idflx  , 'in ')
38*59599516SKenneth E. Jansen          call commu (rmass , ilwork,  1            , 'in ')
39*59599516SKenneth E. Jansen       endif
40*59599516SKenneth E. Jansenc
41*59599516SKenneth E. Jansenc  take care of periodic boundary conditions
42*59599516SKenneth E. Jansenc  but not on surface tension terms in qres(:,10-12)
43*59599516SKenneth E. Jansenc  that are used to compute normal vector
44*59599516SKenneth E. Jansenc
45*59599516SKenneth E. Jansen        idflow = (nflow-1)*nsd
46*59599516SKenneth E. Jansen        do j= 1,nshg
47*59599516SKenneth E. Jansen          if ((btest(iBC(j),10))) then
48*59599516SKenneth E. Jansen            i = iper(j)
49*59599516SKenneth E. Jansen            rmass(i) = rmass(i) + rmass(j)
50*59599516SKenneth E. Jansenc            qres(i,:) = qres(i,:) + qres(j,:)
51*59599516SKenneth E. Jansen            qres(i,1:idflow) = qres(i,1:idflow) + qres(j,1:idflow)
52*59599516SKenneth E. Jansen          endif
53*59599516SKenneth E. Jansen        enddo
54*59599516SKenneth E. Jansen
55*59599516SKenneth E. Jansen        do j= 1,nshg
56*59599516SKenneth E. Jansen          if ((btest(iBC(j),10))) then
57*59599516SKenneth E. Jansen            i = iper(j)
58*59599516SKenneth E. Jansen            rmass(j) = rmass(i)
59*59599516SKenneth E. Jansenc            qres(j,:) = qres(i,:)
60*59599516SKenneth E. Jansen            qres(j,1:idflow) = qres(i,1:idflow)
61*59599516SKenneth E. Jansen          endif
62*59599516SKenneth E. Jansen        enddo
63*59599516SKenneth E. Jansenc
64*59599516SKenneth E. Jansenc.... invert the diagonal mass matrix and find q
65*59599516SKenneth E. Jansenc
66*59599516SKenneth E. Jansen        rmass = one/rmass
67*59599516SKenneth E. Jansen
68*59599516SKenneth E. Jansen       do i=1, idflx
69*59599516SKenneth E. Jansen          qres(:,i) = rmass*qres(:,i)
70*59599516SKenneth E. Jansen       enddo
71*59599516SKenneth E. Jansen       if (isurf .eq. 1) then
72*59599516SKenneth E. Jansen          idflow=(nflow-1)*nsd
73*59599516SKenneth E. Jansenc
74*59599516SKenneth E. Jansenc.... calculation of the unit normal vector
75*59599516SKenneth E. Jansenc
76*59599516SKenneth E. Jansen           tmp =  sqrt(qres(:,idflow+1)**2
77*59599516SKenneth E. Jansen     &               + qres(:,idflow+2)**2
78*59599516SKenneth E. Jansen     &               + qres(:,idflow+3)**2)
79*59599516SKenneth E. Jansen           do i = 1, nshg
80*59599516SKenneth E. Jansen             if (tmp(i) .lt. 0.0001) tmp(i) = 0.0001
81*59599516SKenneth E. Jansen           end do
82*59599516SKenneth E. Jansen           tmp = one/tmp
83*59599516SKenneth E. Jansen
84*59599516SKenneth E. Jansen          do i=1, nsd
85*59599516SKenneth E. Jansen             qres(:,idflow+i) = tmp*qres(:,idflow+i)
86*59599516SKenneth E. Jansen          enddo
87*59599516SKenneth E. Jansen       endif
88*59599516SKenneth E. Jansen
89*59599516SKenneth E. Jansen       if(numpe > 1) then
90*59599516SKenneth E. Jansen          call commu (qres, ilwork, idflx, 'out')
91*59599516SKenneth E. Jansen       endif
92*59599516SKenneth E. Jansen
93*59599516SKenneth E. Jansen       if(iabc==1) then         !are there any axisym bc's
94*59599516SKenneth E. Jansenc
95*59599516SKenneth E. Jansenc       slave has masters value, for abc we need to rotate it
96*59599516SKenneth E. Jansenc
97*59599516SKenneth E. Jansen          do i=1,idflx/nsd
98*59599516SKenneth E. Jansen             do j=1,2
99*59599516SKenneth E. Jansen                istrt=j+(i-1)*(nflow-1)
100*59599516SKenneth E. Jansen                uv(:,j)=qres(:,istrt)
101*59599516SKenneth E. Jansen             enddo
102*59599516SKenneth E. Jansen             call rotabc(uv, iBC, 'out')
103*59599516SKenneth E. Jansen             do j=1,2
104*59599516SKenneth E. Jansen                istrt=j+(i-1)*(nflow-1)
105*59599516SKenneth E. Jansen                qres(:,istrt)=uv(:,j)
106*59599516SKenneth E. Jansen             enddo
107*59599516SKenneth E. Jansen          enddo
108*59599516SKenneth E. Jansen       endif
109*59599516SKenneth E. Jansen
110*59599516SKenneth E. Jansenc
111*59599516SKenneth E. Jansenc.... return
112*59599516SKenneth E. Jansenc
113*59599516SKenneth E. Jansen        return
114*59599516SKenneth E. Jansen        end
115*59599516SKenneth E. Jansen
116*59599516SKenneth E. Jansen
117*59599516SKenneth E. Jansen      subroutine qpbcSclr( rmass, qres, iBC, iper, ilwork )
118*59599516SKenneth E. Jansenc---------------------------------------------------------------------
119*59599516SKenneth E. Jansenc
120*59599516SKenneth E. Jansenc This routine satisfies the periodic boundary conditions
121*59599516SKenneth E. Jansenc on the diffusive flux residual and mass matrix
122*59599516SKenneth E. Jansenc
123*59599516SKenneth E. Jansenc input:
124*59599516SKenneth E. Jansenc     rmass   (nshg)              : mass matrix
125*59599516SKenneth E. Jansenc     qres    (nshg, nsd)         : diffusive flux vector
126*59599516SKenneth E. Jansenc
127*59599516SKenneth E. Jansenc output: modified qres and rmass
128*59599516SKenneth E. Jansenc---------------------------------------------------------------------
129*59599516SKenneth E. Jansen      include "common.h"
130*59599516SKenneth E. Jansen
131*59599516SKenneth E. Jansen      dimension rmass(nshg), qres(nshg,nsd),
132*59599516SKenneth E. Jansen     &          iBC(nshg), iper(nshg)
133*59599516SKenneth E. Jansen
134*59599516SKenneth E. Jansen      if(iabc==1) !are there any axisym bc's
135*59599516SKenneth E. Jansen     &         call rotabc(qres, iBC,  'in ')
136*59599516SKenneth E. Jansen
137*59599516SKenneth E. Jansenc
138*59599516SKenneth E. Jansenc.... compute qi for node A, i.e., qres <-- qres/rmass
139*59599516SKenneth E. Jansenc
140*59599516SKenneth E. Jansen       if (numpe > 1) then
141*59599516SKenneth E. Jansen          call commu (qres  , ilwork, nsd  , 'in ')
142*59599516SKenneth E. Jansen          call commu (rmass , ilwork,  1   , 'in ')
143*59599516SKenneth E. Jansen       endif
144*59599516SKenneth E. Jansen
145*59599516SKenneth E. Jansenc
146*59599516SKenneth E. Jansenc  take care of periodic boundary conditions
147*59599516SKenneth E. Jansenc
148*59599516SKenneth E. Jansen        do j= 1,nshg
149*59599516SKenneth E. Jansen          if (btest(iBC(j),10)) then
150*59599516SKenneth E. Jansen            i = iper(j)
151*59599516SKenneth E. Jansen            rmass(i) = rmass(i) + rmass(j)
152*59599516SKenneth E. Jansen            qres(i,:) = qres(i,:) + qres(j,:)
153*59599516SKenneth E. Jansen          endif
154*59599516SKenneth E. Jansen        enddo
155*59599516SKenneth E. Jansen
156*59599516SKenneth E. Jansen        do j= 1,nshg
157*59599516SKenneth E. Jansen          if (btest(iBC(j),10)) then
158*59599516SKenneth E. Jansen            i = iper(j)
159*59599516SKenneth E. Jansen            rmass(j) = rmass(i)
160*59599516SKenneth E. Jansen            qres(j,:) = qres(i,:)
161*59599516SKenneth E. Jansen          endif
162*59599516SKenneth E. Jansen        enddo
163*59599516SKenneth E. Jansenc
164*59599516SKenneth E. Jansenc.... invert the diagonal mass matrix and find q
165*59599516SKenneth E. Jansenc
166*59599516SKenneth E. Jansen        rmass = one/rmass
167*59599516SKenneth E. Jansen
168*59599516SKenneth E. Jansen       do i=1, nsd
169*59599516SKenneth E. Jansen          qres(:,i) = rmass*qres(:,i)
170*59599516SKenneth E. Jansen       enddo
171*59599516SKenneth E. Jansen
172*59599516SKenneth E. Jansen       if(numpe > 1) then
173*59599516SKenneth E. Jansen          call commu (qres, ilwork, nsd, 'out')
174*59599516SKenneth E. Jansen       endif
175*59599516SKenneth E. Jansen
176*59599516SKenneth E. Jansen      if(iabc==1) !are there any axisym bc's
177*59599516SKenneth E. Jansen     &         call rotabc(qres, iBC, 'out')
178*59599516SKenneth E. Jansen
179*59599516SKenneth E. Jansenc
180*59599516SKenneth E. Jansenc.... return
181*59599516SKenneth E. Jansenc
182*59599516SKenneth E. Jansen        return
183*59599516SKenneth E. Jansen        end
184*59599516SKenneth E. Jansen
185