OpenRadioss 2025.1.11
OpenRadioss project
Loading...
Searching...
No Matches
bemsolvp.F File Reference
#include "implicit_f.inc"
#include "units_c.inc"
#include "task_c.inc"

Go to the source code of this file.

Functions/Subroutines

subroutine bemsolvp (elem, norm, x, nnrp, cr_loc_glob, cr_glob_loc, nncp, cc_loc_glob, cc_glob_loc, nceip, ccei_loc_glob, nep, ce_loc_glob, ce_glob_loc, nreip, crei_loc_glob, nno, nel, hbem, gbem, q, phi, iform, l_assemb, nprow, npcol, nbloc, ipiv, nerp, ilvout, phi_inf)
subroutine glbemp (x1, y1, z1, x2, y2, z2, x3, y3, z3, x0, y0, z0, telnor, hbem, gbem, tbemtg, x, n1, n2, n3, nnrp, nel, jac, cc_glob_loc, ce_glob_loc, nc1, nc2, nc3, offnr, offnc)
subroutine glsinfp (n1, n2, n3, jac, s, phi_inf, nl1, nl2, nl3, offnr)

Function/Subroutine Documentation

◆ bemsolvp()

subroutine bemsolvp ( integer, dimension(3,*) elem,
norm,
x,
integer nnrp,
integer, dimension(*) cr_loc_glob,
integer, dimension(*) cr_glob_loc,
integer nncp,
integer, dimension(*) cc_loc_glob,
integer, dimension(*) cc_glob_loc,
integer nceip,
integer, dimension(*) ccei_loc_glob,
integer nep,
integer, dimension(*) ce_loc_glob,
integer, dimension(*) ce_glob_loc,
integer nreip,
integer, dimension(*) crei_loc_glob,
integer nno,
integer nel,
hbem,
gbem,
q,
phi,
integer iform,
logical l_assemb,
integer nprow,
integer npcol,
integer nbloc,
integer, dimension(*) ipiv,
integer nerp,
integer ilvout,
phi_inf )

Definition at line 35 of file bemsolvp.F.

43C-----------------------------------------------
44C I m p l i c i t T y p e s
45C-----------------------------------------------
46#include "implicit_f.inc"
47C-----------------------------------------------
48C C o m m o n B l o c k s
49C-----------------------------------------------
50#include "units_c.inc"
51#include "task_c.inc"
52C-----------------------------------------------
53C D u m m y A r g u m e n t s
54C-----------------------------------------------
55 INTEGER ELEM(3,*), NNRP, CR_LOC_GLOB(*), CR_GLOB_LOC(*), NNCP,
56 . CC_GLOB_LOC(*), CC_LOC_GLOB(*), NCEIP, CCEI_LOC_GLOB(*),
57 . NEP, CE_LOC_GLOB(*), CE_GLOB_LOC(*), NNO, NEL, IFORM,
58 . NREIP, CREI_LOC_GLOB(*), NPROW, NPCOL, NBLOC, IPIV(*),
59 . NERP, ILVOUT
61 . norm(3,*), x(3,*), hbem(nnrp,*), gbem(nnrp,*), q(*),
62 . phi(*), phi_inf(*)
63 LOGICAL L_ASSEMB
64C-----------------------------------------------
65C L o c a l V a r i a b l e s
66C-----------------------------------------------
67 INTEGER IEL, N1, N2, N3, NL1, NL2, NL3, I, J, JN, IR, IC, JJ,
68 . ICTXT, DESC_H(9), DESC_G(9), DESC_Q(9), DESC_P(9),
69 . INFO, IBID, NC1, NC2, NC3, OFFNR(3), OFFNC(3)
71 . x1, y1, z1, x2, y2, z2, x3, y3, z3, x0, y0, z0, d2,
72 . nrx, nry, nrz, area2, rval(3), lsum(nno), lsumt(nno)
73#if defined(MPI) && defined(MYREAL8) && !defined(WITHOUT_LINALG)
74C---------------------------------------------------
75C 0- Initialisations pour utilisation SCALAPACK
76C---------------------------------------------------
77C Initialisation de la process grid
78 CALL sl_init(ictxt, nprow, npcol)
79C Descripteur pour HBEM
80 desc_h(1)=1
81 desc_h(2)=ictxt
82 desc_h(3)=nno
83 desc_h(4)=nno
84 desc_h(5)=nbloc
85 desc_h(6)=nbloc
86 desc_h(7)=0
87 desc_h(8)=0
88 desc_h(9)=max(1,nnrp)
89C Descripteur pour GBEM
90 desc_g(1)=1
91 desc_g(2)=ictxt
92 desc_g(3)=nno
93 desc_g(4)=nel
94 desc_g(5)=nbloc
95 desc_g(6)=nbloc
96 desc_g(7)=0
97 desc_g(8)=0
98 desc_g(9)=max(1,nnrp)
99C Descripteur pour Q
100 desc_q(1)=1
101 desc_q(2)=ictxt
102 desc_q(3)=nel
103 desc_q(4)=1
104 desc_q(5)=nbloc
105 desc_q(6)=nbloc
106 desc_q(7)=0
107 desc_q(8)=0
108 desc_q(9)=max(1,nerp)
109C Descripteur pour PHI
110 desc_p(1)=1
111 desc_p(2)=ictxt
112 desc_p(3)=nno
113 desc_p(4)=1
114 desc_p(5)=nbloc
115 desc_p(6)=nbloc
116 desc_p(7)=0
117 desc_p(8)=0
118 desc_p(9)=max(1,nnrp)
119C---------------------------------------------------
120C 1- Assemblage
121C---------------------------------------------------
122 IF (l_assemb) THEN
123 IF (ilvout>=1.AND.ispmd==0) WRITE(istdo,'(A)') ' ** BEMSOLV : ASSEMBLY OF INTEGRAL OPERATORS'
124 IF (iform==1) THEN
125C Assemblage de H
126 DO i=1,nnrp
127 DO j=1,nncp
128 hbem(i,j)=zero
129 ENDDO
130 DO j=1,nep
131 gbem(i,j)=zero
132 ENDDO
133 ENDDO
134 DO i=1,nceip
135 IF (ilvout==2)
136 . CALL progcondp_c(i,nceip+nep, ispmd+1, ibid)
137 iel=ccei_loc_glob(i)
138 n1=elem(1,iel)
139 n2=elem(2,iel)
140 n3=elem(3,iel)
141 nl1=cc_glob_loc(n1)
142 nl2=cc_glob_loc(n2)
143 nl3=cc_glob_loc(n3)
144 x1=x(1,n1)
145 x2=x(1,n2)
146 x3=x(1,n3)
147 y1=x(2,n1)
148 y2=x(2,n2)
149 y3=x(2,n3)
150 z1=x(3,n1)
151 z2=x(3,n2)
152 z3=x(3,n3)
153 x0=third*(x1+x2+x3)
154 y0=third*(y1+y2+y3)
155 z0=third*(z1+z2+z3)
156 d2=min((x0-x1)**2+(y0-y1)**2+(z0-z1)**2,
157 . (x0-x2)**2+(y0-y2)**2+(z0-z2)**2,
158 . (x0-x3)**2+(y0-y3)**2+(z0-z3)**2)
159 nrx=norm(1,iel)
160 nry=norm(2,iel)
161 nrz=norm(3,iel)
162 area2=sqrt(nrx**2+nry**2+nrz**2)
163 DO j=1,nnrp
164 jn=cr_loc_glob(j)
165 CALL inthtg(jn, x1 , y1, z1, x2,
166 . y2, z2, x3, y3, z3,
167 . x0, y0, z0, d2, area2,
168 . nrx, nry, nrz, rval, n1,
169 . n2, n3, x )
170 IF (nl1/=0) hbem(j,nl1)=hbem(j,nl1)+rval(1)
171 IF (nl2/=0) hbem(j,nl2)=hbem(j,nl2)+rval(2)
172 IF (nl3/=0) hbem(j,nl3)=hbem(j,nl3)+rval(3)
173 ENDDO
174 ENDDO
175C Assemblage de G
176 DO i=1,nep
177 IF (ilvout==2)
178 . CALL progcondp_c(nceip+i,nceip+nep, ispmd+1, ibid)
179 iel=ce_loc_glob(i)
180 n1=elem(1,iel)
181 n2=elem(2,iel)
182 n3=elem(3,iel)
183 nl1=cc_glob_loc(n1)
184 nl2=cc_glob_loc(n2)
185 nl3=cc_glob_loc(n3)
186 x1=x(1,n1)
187 x2=x(1,n2)
188 x3=x(1,n3)
189 y1=x(2,n1)
190 y2=x(2,n2)
191 y3=x(2,n3)
192 z1=x(3,n1)
193 z2=x(3,n2)
194 z3=x(3,n3)
195 x0=third*(x1+x2+x3)
196 y0=third*(y1+y2+y3)
197 z0=third*(z1+z2+z3)
198 d2=min((x0-x1)**2+(y0-y1)**2+(z0-z1)**2,
199 . (x0-x2)**2+(y0-y2)**2+(z0-z2)**2,
200 . (x0-x3)**2+(y0-y3)**2+(z0-z3)**2)
201 nrx=norm(1,iel)
202 nry=norm(2,iel)
203 nrz=norm(3,iel)
204 area2=sqrt(nrx**2+nry**2+nrz**2)
205 DO j=1,nnrp
206 jn=cr_loc_glob(j)
207 CALL intgtg(jn, x1 , y1, z1, x2,
208 . y2, z2, x3, y3, z3,
209 . x0, y0, z0, d2, area2,
210 . nrx, nry, nrz, rval, n1,
211 . n2, n3, x )
212 gbem(j,i)=gbem(j,i)+rval(1)
213 ENDDO
214 ENDDO
215C Somme des lignes pour le terme diagonal de H
216 DO i=1,nno
217 lsum(i)=zero
218 ENDDO
219 DO i=1,nncp
220 DO j=1,nnrp
221 jj=cr_loc_glob(j)
222 lsum(jj)=lsum(jj)+hbem(j,i)
223 ENDDO
224 ENDDO
225 CALL spmd_fl_sum(lsum, nno, lsumt)
226 DO i=1,nno
227 ir=cr_glob_loc(i)
228 ic=cc_glob_loc(i)
229 IF (ir>0.AND.ic>0) hbem(ir,ic)=-lsumt(i)
230 ENDDO
231C
232 ELSEIF (iform==2) THEN
233 DO i=1,nnrp
234 DO j=1,nncp
235 hbem(i,j)=zero
236 ENDDO
237 DO j=1,nep
238 gbem(i,j)=zero
239 ENDDO
240 ENDDO
241 DO i=1,nreip
242 IF (ilvout==2)
243 . CALL progcondp_c(i,nreip, ispmd+1, ibid)
244 iel=crei_loc_glob(i)
245 n1=elem(1,iel)
246 n2=elem(2,iel)
247 n3=elem(3,iel)
248C
249 nl1=cr_glob_loc(n1)
250 nl2=cr_glob_loc(n2)
251 nl3=cr_glob_loc(n3)
252 offnr(1)=min(1,nl1)
253 offnr(2)=min(1,nl2)
254 offnr(3)=min(1,nl3)
255 nl1=max(1,nl1)
256 nl2=max(1,nl2)
257 nl3=max(1,nl3)
258C
259 nc1=cc_glob_loc(n1)
260 nc2=cc_glob_loc(n2)
261 nc3=cc_glob_loc(n3)
262 offnc(1)=min(1,nc1)
263 offnc(2)=min(1,nc2)
264 offnc(3)=min(1,nc3)
265 nc1=max(1,nc1)
266 nc2=max(1,nc2)
267 nc3=max(1,nc3)
268C
269 x1=x(1,n1)
270 x2=x(1,n2)
271 x3=x(1,n3)
272 y1=x(2,n1)
273 y2=x(2,n2)
274 y3=x(2,n3)
275 z1=x(3,n1)
276 z2=x(3,n2)
277 z3=x(3,n3)
278 x0=third*(x1+x2+x3)
279 y0=third*(y1+y2+y3)
280 z0=third*(z1+z2+z3)
281 nrx=norm(1,iel)
282 nry=norm(2,iel)
283 nrz=norm(3,iel)
284 area2=sqrt(nrx**2+nry**2+nrz**2)
285 CALL glbemp(x1, y1, z1, x2, y2,
286 . z2, x3, y3, z3, x0,
287 . y0, z0, norm, hbem, gbem,
288 . elem, x, nl1, nl2, nl3,
289 . nnrp, nel, area2, cc_glob_loc,
290 . ce_glob_loc, nc1, nc2, nc3 ,
291 . offnr, offnc)
292 ENDDO
293 ENDIF
294C Somme des lignes de G pour remplacer la derniere colonne de H
295 DO i=1,nno
296 lsum(i)=zero
297 ENDDO
298 DO i=1,nep
299 DO j=1,nnrp
300 jj=cr_loc_glob(j)
301 lsum(jj)=lsum(jj)+gbem(j,i)
302 ENDDO
303 ENDDO
304 CALL spmd_fl_sum(lsum, nno, lsumt)
305 ic=cc_glob_loc(nno)
306 IF (ic>0) THEN
307 DO i=1,nno
308 ir=cr_glob_loc(i)
309 IF (ir>0) hbem(ir,ic)=-lsumt(i)
310 ENDDO
311 ENDIF
312C Factorisation de HBEM
313 CALL pdgetrf(nno, nno, hbem, 1, 1,
314 . desc_h, ipiv, info)
315 ENDIF
316C---------------------------------------------------
317C 2- Resolution
318C---------------------------------------------------
319 IF (ilvout>=1.AND.ispmd==0) THEN
320 WRITE(istdo,*)
321 WRITE(istdo,'(A)') ' ** BEMSOLV : PARALLEL LINEAR SYSTEM SOLUTION'
322 ENDIF
323 CALL pdgemv('N', nno, nel, one, gbem,
324 . 1, 1, desc_g, q, 1,
325 . 1, desc_q, 1, zero, phi,
326 . 1, 1, desc_p, 1)
327C
328 IF (iform==1) THEN
329 ic=cc_glob_loc(1)
330 IF (ic>0) THEN
331 DO i=1,nnrp
332 ir=cr_loc_glob(i)
333 phi(i)=phi(i)-phi_inf(ir)
334 ENDDO
335 ENDIF
336 ELSEIF (iform==2) THEN
337 ic=cc_glob_loc(1)
338 IF (ic>0) THEN
339 DO i=1,nreip
340 iel=crei_loc_glob(i)
341 n1=elem(1,iel)
342 n2=elem(2,iel)
343 n3=elem(3,iel)
344C
345 nl1=cr_glob_loc(n1)
346 nl2=cr_glob_loc(n2)
347 nl3=cr_glob_loc(n3)
348 offnr(1)=min(1,nl1)
349 offnr(2)=min(1,nl2)
350 offnr(3)=min(1,nl3)
351 nl1=max(1,nl1)
352 nl2=max(1,nl2)
353 nl3=max(1,nl3)
354C
355 nrx=norm(1,iel)
356 nry=norm(2,iel)
357 nrz=norm(3,iel)
358 area2=sqrt(nrx**2+nry**2+nrz**2)
359 CALL glsinfp(n1, n2, n3, area2, phi,
360 . phi_inf, nl1, nl2, nl3 , offnr)
361 ENDDO
362 ENDIF
363 ENDIF
364C
365 CALL pdgetrs('N', nno, 1, hbem, 1,
366 . 1, desc_h, ipiv, phi, 1,
367 . 1, desc_p, info)
368C
369 CALL blacs_gridexit(ictxt)
370C
371#endif
372 RETURN
subroutine sl_init(ictxt, nprow, npcol)
Definition SL_init.f:2
subroutine glsinfp(n1, n2, n3, jac, s, phi_inf, nl1, nl2, nl3, offnr)
Definition bemsolvp.F:548
subroutine glbemp(x1, y1, z1, x2, y2, z2, x3, y3, z3, x0, y0, z0, telnor, hbem, gbem, tbemtg, x, n1, n2, n3, nnrp, nel, jac, cc_glob_loc, ce_glob_loc, nc1, nc2, nc3, offnr, offnc)
Definition bemsolvp.F:388
#define my_real
Definition cppsort.cpp:32
norm(diag(diag(diag(inv(mat))) -id.SOL), 2) % destroy mumps instance id.JOB
subroutine intgtg(x1, y1, z1, x2, y2, z2, x3, y3, z3, xp, yp, zp, d2, jac, xs, ys, zs, rval)
Definition inte_tg.F:145
subroutine inthtg(x1, y1, z1, x2, y2, z2, x3, y3, z3, xp, yp, zp, nrx, nry, nrz, d2, jac, xs, ys, zs, rval)
Definition inte_tg.F:33
#define min(a, b)
Definition macros.h:20
#define max(a, b)
Definition macros.h:21
subroutine pdgetrs(trans, n, nrhs, a, ia, ja, desca, ipiv, b, ib, jb, descb, info)
Definition mpi.f:1171
subroutine pdgetrf(m, n, a, ia, ja, desca, ipiv, info)
Definition mpi.f:914
subroutine blacs_gridexit(cntxt)
Definition mpi.f:762
void progcondp_c(int *icur, int *imax, int *iproc, int *id)
Definition progcond_c.c:61
subroutine spmd_fl_sum(lsum, len, lsumt)
Definition spmd_fl_sum.F:37

◆ glbemp()

subroutine glbemp ( x1,
y1,
z1,
x2,
y2,
z2,
x3,
y3,
z3,
x0,
y0,
z0,
telnor,
hbem,
gbem,
integer, dimension(3,*) tbemtg,
x,
integer n1,
integer n2,
integer n3,
integer nnrp,
integer nel,
jac,
integer, dimension(*) cc_glob_loc,
integer, dimension(*) ce_glob_loc,
integer nc1,
integer nc2,
integer nc3,
integer, dimension(*) offnr,
integer, dimension(*) offnc )

Definition at line 381 of file bemsolvp.F.

388C-----------------------------------------------
389C I m p l i c i t T y p e s
390C-----------------------------------------------
391#include "implicit_f.inc"
392C-----------------------------------------------
393C D u m m y A r g u m e n t s
394C-----------------------------------------------
395 INTEGER TBEMTG(3,*), N1, N2, N3, NEL, NNRP, CC_GLOB_LOC(*),
396 . CE_GLOB_LOC(*), NC1, NC2, NC3, OFFNR(*), OFFNC(*)
397 my_real
398 . x1, y1, z1, x2, y2, z2, x3, y3, z3, x0, y0, z0,
399 . telnor(3,*), x(3,*), jac, hbem(nnrp,*), gbem(nnrp,*)
400C-----------------------------------------------
401C L o c a l V a r i a b l e s
402C-----------------------------------------------
403 INTEGER NPG, IAD, JBID, NBID, IDEG, IAD2, IP, IEL, NN1, NN2, NN3,
404 . NL1, NL2, NL3, IELL, OFFNR2(3), OFFEL
405 my_real
406 . pg(50), wpg(25), w, ksip, etap, val1, val2, val3, xp,
407 . yp, zp, cp, xx1, yy1, zz1, xx2, yy2, zz2, xx3, yy3, zz3,
408 . xx4, yy4, zz4, xx5, yy5, zz5, xx6, yy6, zz6, xx0, yy0,
409 . zz0, nrx, nry, nrz, d2, rvl(6), rvlh(3), rvlg
410 .
411C
412 DATA pg /.33333333,.33333333,
413 . .33333333,.33333333,
414 . .60000000,.20000000,
415 . .20000000,.60000000,
416 . .20000000,.20000000,
417 . .33333333,.33333333,
418 . .79742699,.10128651,
419 . .10128651,.79742699,
420 . .10128651,.10128651,
421 . .05971587,.47014206,
422 . .47014206,.05971587,
423 . .47014206,.47014206,
424 . .06513010,.06513010,
425 . .86973979,.06513010,
426 . .06513010,.86973979,
427 . .31286550,.04869031,
428 . .63844419,.31286550,
429 . .04869031,.63844419,
430 . .63844419,.04869031,
431 . .31286550,.63844419,
432 . .04869031,.31286550,
433 . .26034597,.26034597,
434 . .47930807,.26034597,
435 . .26034597,.47930807,
436 . .33333333,.33333333/
437 DATA wpg /1.00000000,
438 . -.56250000,.52083333,
439 . .52083333,.52083333,
440 . .22500000,.12593918,
441 . .12593918,.12593918,
442 . .13239415,.13239415,
443 . .13239415,
444 . .05334724,.05334724,
445 . .05334724,.07711376,
446 . .07711376,.07711376,
447 . .07711376,.07711376,
448 . .07711376,.17561526,
449 . .17561526,.17561526,
450 . -.14957004/
451C
452 npg=13
453 iad=13
454C
455 iad2=2*(iad-1)+1
456 DO ip=1,npg
457 w=wpg(iad)
458 ksip=pg(iad2)
459 etap=pg(iad2+1)
460 iad=iad+1
461 iad2=iad2+2
462 val1=one-ksip-etap
463 val2=ksip
464 val3=etap
465 xp=val1*x1+val2*x2+val3*x3
466 yp=val1*y1+val2*y2+val3*y3
467 zp=val1*z1+val2*z2+val3*z3
468 cp=zero
469 DO iel=1,nel
470 iell=ce_glob_loc(iel)
471 offel=min(1,iell)
472 iell=max(1,iell)
473C
474 nn1=tbemtg(1,iel)
475 nn2=tbemtg(2,iel)
476 nn3=tbemtg(3,iel)
477 nl1=cc_glob_loc(nn1)
478 nl2=cc_glob_loc(nn2)
479 nl3=cc_glob_loc(nn3)
480 offnr2(1)=min(1,nl1)
481 offnr2(2)=min(1,nl2)
482 offnr2(3)=min(1,nl3)
483 nl1=max(1,nl1)
484 nl2=max(1,nl2)
485 nl3=max(1,nl3)
486C
487 xx1=x(1,nn1)
488 yy1=x(2,nn1)
489 zz1=x(3,nn1)
490 xx2=x(1,nn2)
491 yy2=x(2,nn2)
492 zz2=x(3,nn2)
493 xx3=x(1,nn3)
494 yy3=x(2,nn3)
495 zz3=x(3,nn3)
496 nrx=telnor(1,iel)
497 nry=telnor(2,iel)
498 nrz=telnor(3,iel)
499 CALL intanl(xx1 , yy1 , zz1, xx2, yy2,
500 . zz2 , xx3 , yy3, zz3, xp ,
501 . yp , zp , nrx, nry, nrz,
502 . rvlh, rvlg)
503C Matrice H
504 hbem(n1,nl1)=hbem(n1,nl1)
505 . +offnr(1)*offnr2(1)*w*val1*rvlh(1)*jac
506 hbem(n1,nl2)=hbem(n1,nl2)
507 . +offnr(1)*offnr2(2)*w*val1*rvlh(2)*jac
508 hbem(n1,nl3)=hbem(n1,nl3)
509 . +offnr(1)*offnr2(3)*w*val1*rvlh(3)*jac
510 hbem(n2,nl1)=hbem(n2,nl1)
511 . +offnr(2)*offnr2(1)*w*val2*rvlh(1)*jac
512 hbem(n2,nl2)=hbem(n2,nl2)
513 . +offnr(2)*offnr2(2)*w*val2*rvlh(2)*jac
514 hbem(n2,nl3)=hbem(n2,nl3)
515 . +offnr(2)*offnr2(3)*w*val2*rvlh(3)*jac
516 hbem(n3,nl1)=hbem(n3,nl1)
517 . +offnr(3)*offnr2(1)*w*val3*rvlh(1)*jac
518 hbem(n3,nl2)=hbem(n3,nl2)
519 . +offnr(3)*offnr2(2)*w*val3*rvlh(2)*jac
520 hbem(n3,nl3)=hbem(n3,nl3)
521 . +offnr(3)*offnr2(3)*w*val3*rvlh(3)*jac
522 cp=cp-rvlh(1)-rvlh(2)-rvlh(3)
523C Matrice G
524 gbem(n1,iell)=gbem(n1,iell)+offnr(1)*offel*w*val1*rvlg*jac
525 gbem(n2,iell)=gbem(n2,iell)+offnr(2)*offel*w*val2*rvlg*jac
526 gbem(n3,iell)=gbem(n3,iell)+offnr(3)*offel*w*val3*rvlg*jac
527 ENDDO
528 hbem(n1,nc1)=hbem(n1,nc1)+offnr(1)*offnc(1)*w*cp*val1*val1*jac
529 hbem(n1,nc2)=hbem(n1,nc2)+offnr(1)*offnc(2)*w*cp*val1*val2*jac
530 hbem(n1,nc3)=hbem(n1,nc3)+offnr(1)*offnc(3)*w*cp*val1*val3*jac
531 hbem(n2,nc1)=hbem(n2,nc1)+offnr(2)*offnc(1)*w*cp*val2*val1*jac
532 hbem(n2,nc2)=hbem(n2,nc2)+offnr(2)*offnc(2)*w*cp*val2*val2*jac
533 hbem(n2,nc3)=hbem(n2,nc3)+offnr(2)*offnc(3)*w*cp*val2*val3*jac
534 hbem(n3,nc1)=hbem(n3,nc1)+offnr(3)*offnc(1)*w*cp*val3*val1*jac
535 hbem(n3,nc2)=hbem(n3,nc2)+offnr(3)*offnc(2)*w*cp*val3*val2*jac
536 hbem(n3,nc3)=hbem(n3,nc3)+offnr(3)*offnc(3)*w*cp*val3*val3*jac
537 ENDDO
538C
539 RETURN
subroutine intanl(x1, y1, z1, x2, y2, z2, x3, y3, z3, xs, ys, zs, nrx, nry, nrz, rvlh, rvlg)
Definition bemsolv.F:733

◆ glsinfp()

subroutine glsinfp ( integer n1,
integer n2,
integer n3,
jac,
s,
phi_inf,
integer nl1,
integer nl2,
integer nl3,
integer, dimension(*) offnr )

Definition at line 546 of file bemsolvp.F.

548C-----------------------------------------------
549C I m p l i c i t T y p e s
550C-----------------------------------------------
551#include "implicit_f.inc"
552C-----------------------------------------------
553C D u m m y A r g u m e n t s
554C-----------------------------------------------
555 INTEGER N1, N2, N3, NL1, NL2, NL3, OFFNR(*)
556 my_real
557 . jac, s(*), phi_inf(*)
558C-----------------------------------------------
559C L o c a l V a r i a b l e s
560C-----------------------------------------------
561 INTEGER NPG, IAD, IAD2, IP
562 my_real
563 . pg(50), wpg(25), w, ksip, etap, val1, val2, val3,
564 . phi_infp
565C
566 DATA pg /.33333333,.33333333,
567 . .33333333,.33333333,
568 . .60000000,.20000000,
569 . .20000000,.60000000,
570 . .20000000,.20000000,
571 . .33333333,.33333333,
572 . .79742699,.10128651,
573 . .10128651,.79742699,
574 . .10128651,.10128651,
575 . .05971587,.47014206,
576 . .47014206,.05971587,
577 . .47014206,.47014206,
578 . .06513010,.06513010,
579 . .86973979,.06513010,
580 . .06513010,.86973979,
581 . .31286550,.04869031,
582 . .63844419,.31286550,
583 . .04869031,.63844419,
584 . .63844419,.04869031,
585 . .31286550,.63844419,
586 . .04869031,.31286550,
587 . .26034597,.26034597,
588 . .47930807,.26034597,
589 . .26034597,.47930807,
590 . .33333333,.33333333/
591 DATA wpg /1.00000000,
592 . -.56250000,.52083333,
593 . .52083333,.52083333,
594 . .22500000,.12593918,
595 . .12593918,.12593918,
596 . .13239415,.13239415,
597 . .13239415,
598 . .05334724,.05334724,
599 . .05334724,.07711376,
600 . .07711376,.07711376,
601 . .07711376,.07711376,
602 . .07711376,.17561526,
603 . .17561526,.17561526,
604 . -.14957004/
605C
606 npg=13
607 iad=13
608C
609 iad2=2*(iad-1)+1
610 DO ip=1,npg
611 w=wpg(iad)
612 ksip=pg(iad2)
613 etap=pg(iad2+1)
614 iad=iad+1
615 iad2=iad2+2
616 val1=one-ksip-etap
617 val2=ksip
618 val3=etap
619 phi_infp=phi_inf(n1)*val1+phi_inf(n2)*val2+phi_inf(n3)*val3
620 s(nl1)=s(nl1)-offnr(1)*w*val1*phi_infp*jac
621 s(nl2)=s(nl2)-offnr(2)*w*val2*phi_infp*jac
622 s(nl3)=s(nl3)-offnr(3)*w*val3*phi_infp*jac
623* IF (NL1>0) S(NL1)=S(NL1)-W*VAL1*PHI_INFP*JAC
624* IF (NL2>0) S(NL2)=S(NL2)-W*VAL2*PHI_INFP*JAC
625* IF (NL3>0) S(NL3)=S(NL3)-W*VAL3*PHI_INFP*JAC
626 ENDDO
627C
628 RETURN