OpenRadioss 2025.1.11
OpenRadioss project
Loading...
Searching...
No Matches
c3fint_reg_ini.F File Reference
#include "implicit_f.inc"
#include "param_c.inc"
#include "com01_c.inc"
#include "com04_c.inc"
#include "scr03_c.inc"
#include "vectorize.inc"

Go to the source code of this file.

Functions/Subroutines

subroutine c3fint_reg_ini (elbuf_tab, nloc_dmg, area, ixtg, dt_nl, x, xreftg, nft, nel, ng, ipm, bufmat, time, failure)

Function/Subroutine Documentation

◆ c3fint_reg_ini()

subroutine c3fint_reg_ini ( type(elbuf_struct_), dimension(ngroup), target elbuf_tab,
type (nlocal_str_), target nloc_dmg,
intent(in) area,
integer, dimension(nixtg,*) ixtg,
dt_nl,
x,
xreftg,
integer nft,
integer nel,
integer ng,
integer, dimension(npropmi,*) ipm,
bufmat,
time,
logical failure )

Definition at line 33 of file c3fint_reg_ini.F.

37C-----------------------------------------------
38C M o d u l e s
39C-----------------------------------------------
40 USE elbufdef_mod
41 USE message_mod
43 use element_mod , only : nixtg
44C-----------------------------------------------
45C I m p l i c i t T y p e s
46C-----------------------------------------------
47#include "implicit_f.inc"
48C-----------------------------------------------
49C C o m m o n B l o c k s
50C-----------------------------------------------
51#include "param_c.inc"
52#include "com01_c.inc"
53#include "com04_c.inc"
54#include "scr03_c.inc"
55C-----------------------------------------------
56C D u m m y A r g u m e n t s
57C-----------------------------------------------
58 TYPE(ELBUF_STRUCT_), TARGET, DIMENSION(NGROUP) :: ELBUF_TAB
59 TYPE (NLOCAL_STR_) , TARGET :: NLOC_DMG
60 INTEGER IXTG(NIXTG,*),NFT,NEL,NG,IPM(NPROPMI,*)
61 my_real ,DIMENSION(NUMELC+NUMELTG),INTENT(IN) ::
62 . area
64 . x(3,*),xreftg(3,3,*),dt_nl,bufmat(*),time
65 LOGICAL :: FAILURE
66C-----------------------------------------------
67C L o c a l V a r i a b l e s
68C-----------------------------------------------
69 INTEGER :: IMAT,NDOF,L_NLOC,N1,N2,N3,K,I,NDNOD
70 INTEGER, DIMENSION(NEL) :: POS1,POS2,POS3
72 . le_min,len,damp,dens,ntn_unl,ntn_vnl,
73 . ntvar,z01(11,11),wf1(11,11),zn1(12,11),b1,b2,b3,
74 . nth1,nth2,bth1,bth2,k1,k12,k2,sspnl,le_max
75 my_real, DIMENSION(:,:), ALLOCATABLE :: var_reg,vpred
76 my_real,
77 . DIMENSION(:), POINTER :: fnl,unl,vnl,dnl,mnl,thck
78 my_real, DIMENSION(NEL) :: x1,x2,x3,y1,y2,y3,
79 . px1,py1,py2,e1x,e2x,e3x,e1y,e2y,e3y,e1z,e2z,e3z,
80 . x2l,x3l,y3l,z1,z2,z3,surf,x31,y31,z31,offg,
81 . vols,btb11,btb12,btb13,btb22,btb23,btb33
82 TYPE(BUF_NLOC_),POINTER :: BUFNL
83 my_real, DIMENSION(:,:), POINTER ::
84 . massth,fnlth,vnlth,unlth
85 ! Position of integration points in the thickness
86 DATA z01/
87 1 0. ,0. ,0. ,0. ,0. ,
88 1 0. ,0. ,0. ,0. ,0. ,0. ,
89 2 -.5 ,0.5 ,0. ,0. ,0. ,
90 2 0. ,0. ,0. ,0. ,0. ,0. ,
91 3 -.5 ,0. ,0.5 ,0. ,0. ,
92 3 0. ,0. ,0. ,0. ,0. ,0. ,
93 4 -.5 ,-.1666667,0.1666667,0.5 ,0. ,
94 4 0. ,0. ,0. ,0. ,0. ,0. ,
95 5 -.5 ,-.25 ,0. ,0.25 ,0.5 ,
96 5 0. ,0. ,0. ,0. ,0. ,0. ,
97 6 -.5 ,-.3 ,-.1 ,0.1 ,0.3 ,
98 6 0.5 ,0. ,0. ,0. ,0. ,0. ,
99 7 -.5 ,-.3333333,-.1666667,0.0 ,0.1666667,
100 7 0.3333333,0.5 ,0. ,0. ,0. ,0. ,
101 8 -.5 ,-.3571429,-.2142857,-.0714286,0.0714286,
102 8 0.2142857,0.3571429,0.5 ,0. ,0. ,0. ,
103 9 -.5 ,-.375 ,-.25 ,-.125 ,0.0 ,
104 9 0.125 ,0.25 ,0.375 ,0.5 ,0. ,0. ,
105 a -.5 ,-.3888889,-.2777778,-.1666667,-.0555555,
106 a 0.0555555,0.1666667,0.2777778,0.3888889,0.5 ,0. ,
107 b -.5 ,-.4 ,-.3 ,-.2 ,-.1 ,
108 b 0. ,0.1 ,0.2 ,0.3 ,0.4 ,0.5 /
109 ! Weight of integration in the thickness
110 DATA wf1/
111 1 1. ,0. ,0. ,0. ,0. ,
112 1 0. ,0. ,0. ,0. ,0. ,0. ,
113 2 0.5 ,0.5 ,0. ,0. ,0. ,
114 2 0. ,0. ,0. ,0. ,0. ,0. ,
115 3 0.25 ,0.5 ,0.25 ,0. ,0. ,
116 3 0. ,0. ,0. ,0. ,0. ,0. ,
117 4 0.1666667,0.3333333,0.3333333,0.1666667,0. ,
118 4 0. ,0. ,0. ,0. ,0. ,0. ,
119 5 0.125 ,0.25 ,0.25 ,0.25 ,0.125 ,
120 5 0. ,0. ,0. ,0. ,0. ,0. ,
121 6 0.1 ,0.2 ,0.2 ,0.2 ,0.2 ,
122 6 0.1 ,0. ,0. ,0. ,0. ,0. ,
123 7 0.0833333,0.1666667,0.1666667,0.1666667,0.1666667,
124 7 0.1666667,0.0833333,0. ,0. ,0. ,0. ,
125 8 0.0714286,0.1428571,0.1428571,0.1428571,0.1428571,
126 8 0.1428571,0.1428571,0.0714286,0. ,0. ,0. ,
127 9 0.0625 ,0.125 ,0.125 ,0.125 ,0.125 ,
128 9 0.125 ,0.125 ,0.125 ,0.0625 ,0. ,0. ,
129 a 0.0555556,0.1111111,0.1111111,0.1111111,0.1111111,
130 a 0.1111111,0.1111111,0.1111111,0.1111111,0.0555556,0. ,
131 b 0.05 ,0.1 ,0.1 ,0.1 ,0.1 ,
132 b 0.1 ,0.1 ,0.1 ,0.1 ,0.1 ,0.05 /
133 ! Position of nodes in the shell thickness
134 DATA zn1/
135 1 0. ,0. ,0. ,0. ,0. ,0. ,
136 1 0. ,0. ,0. ,0. ,0. ,0. ,
137 2 -.5 ,0.5 ,0. ,0. ,0. ,0. ,
138 2 0. ,0. ,0. ,0. ,0. ,0. ,
139 3 -.5 ,-.25 ,0.25 ,0.5 ,0. ,0. ,
140 3 0. ,0. ,0. ,0. ,0. ,0. ,
141 4 -.5 ,-.3333333,0. ,0.3333333,0.5 ,0. ,
142 4 0. ,0. ,0. ,0. ,0. ,0. ,
143 5 -.5 ,-.375 ,-0.125 ,0.125 ,0.375 ,0.5 ,
144 5 0. ,0. ,0. ,0. ,0. ,0. ,
145 6 -.5 ,-.4 ,-.2 ,0.0 ,0.2 ,0.4 ,
146 6 0.5 ,0. ,0. ,0. ,0. ,0. ,
147 7 -.5 ,-.4166667,-.25 ,-.0833333,0.0833333,0.25 ,
148 7 0.4166667,0.5 ,0. ,0. ,0. ,0. ,
149 8 -.5 ,-.4285715,-.2857143,-.1428572,0.0 ,0.1428572,
150 8 0.2857143,0.4285715,0.5 ,0. ,0. ,0. ,
151 9 -.5 ,-.4375 ,-.3125 ,-.1875 ,-.0625 ,0.0625 ,
152 9 0.1875 ,0.3125 ,0.4375 ,0.5 ,0. ,0. ,
153 a -.5 ,-.4444444,-.3333333,-.2222222,-.1111111,0. ,
154 a 0.1111111,0.2222222,0.3333333,0.4444444,0.5 ,0. ,
155 b -.5 ,-.45 ,-.35 ,-.25 ,-.15 ,-.05 ,
156 b 0.05 ,0.15 ,0.25 ,0.35 ,0.45 ,0.5 /
157C=======================================================================
158 ! Size of the non-local vectors
159 l_nloc = nloc_dmg%L_NLOC
160 ! Pointing the non-local forces vector
161 fnl => nloc_dmg%FNL(1:l_nloc,1)
162 vnl => nloc_dmg%VNL(1:l_nloc)
163 dnl => nloc_dmg%DNL(1:l_nloc)
164 unl => nloc_dmg%UNL(1:l_nloc)
165 mnl => nloc_dmg%MASS(1:l_nloc)
166 ! Number of the material law
167 imat = ixtg(1,1+nft)
168 ! Minimal length in the surface
169 le_min = sqrt((four/sqrt(three))*minval(area(numelc+nft+1:numelc+nft+nel)))
170 ! Number of integration points
171 ndof = elbuf_tab(ng)%BUFLY(1)%NPTT
172 ! Thickness of the shell
173 thck => elbuf_tab(ng)%GBUF%THK(1:nel)
174 ! Global minimal length
175 IF (ndof>1) THEN
176 le_min = min(le_min,minval(thck(1:nel))/ndof)
177 ENDIF
178 ! Non-local internal length
179 len = nloc_dmg%LEN(imat)
180 ! Maximal length of convergence
181 le_max = nloc_dmg%LE_MAX(imat)
182 ! Non-local damping
183 damp = nloc_dmg%DAMP(imat)
184 ! Non-local density
185 dens = nloc_dmg%DENS(imat)
186 ! Non-local sound speed
187 sspnl = nloc_dmg%SSPNL(imat)
188 ! Non-local timestep
189 dt_nl = min(dt_nl,0.5d0*((two*min(le_min,le_max)*sqrt(three*dens))/
190 . (sqrt(twelve*(len**2)+(min(le_min,le_max)**2)))))
191 ! Allocation of the velocities predictor
192 IF (ndof>1) THEN
193 IF (ndof > 2) THEN
194 ALLOCATE(vpred(nel,ndof+1))
195 ndnod = ndof + 1
196 ELSE
197 ALLOCATE(vpred(nel,ndof))
198 ndnod = ndof
199 ENDIF
200 ENDIF
201 ! Variable to regularize
202 IF (.NOT.ALLOCATED(var_reg)) ALLOCATE(var_reg(nel,ndof))
203c
204 ! Computation of kinematical data
205# include "vectorize.inc"
206 DO i = 1,nel
207 ! Coordinates of the nodes of the element
208 IF (nxref == 0) THEN
209 x1(i)=x(1,ixtg(2,nft+i))
210 y1(i)=x(2,ixtg(2,nft+i))
211 z1(i)=x(3,ixtg(2,nft+i))
212 x2(i)=x(1,ixtg(3,nft+i))
213 y2(i)=x(2,ixtg(3,nft+i))
214 z2(i)=x(3,ixtg(3,nft+i))
215 x3(i)=x(1,ixtg(4,nft+i))
216 y3(i)=x(2,ixtg(4,nft+i))
217 z3(i)=x(3,ixtg(4,nft+i))
218 ELSE
219 x1(i)=xreftg(1,1,nft+i)
220 y1(i)=xreftg(1,2,nft+i)
221 z1(i)=xreftg(1,3,nft+i)
222 x2(i)=xreftg(2,1,nft+i)
223 y2(i)=xreftg(2,2,nft+i)
224 z2(i)=xreftg(2,3,nft+i)
225 x3(i)=xreftg(3,1,nft+i)
226 y3(i)=xreftg(3,2,nft+i)
227 z3(i)=xreftg(3,3,nft+i)
228 ENDIF
229C
230 ! Recovering the nodes of the non-local element
231 n1 = nloc_dmg%IDXI(ixtg(2,nft+i))
232 n2 = nloc_dmg%IDXI(ixtg(3,nft+i))
233 n3 = nloc_dmg%IDXI(ixtg(4,nft+i))
234 ! Recovering the positions of the first d.o.fs of each nodes
235 pos1(i) = nloc_dmg%POSI(n1)
236 pos2(i) = nloc_dmg%POSI(n2)
237 pos3(i) = nloc_dmg%POSI(n3)
238 ENDDO
239c
240 ! Loop over Gauss points in the thickness
241 DO k = 1,ndof
242 ! Loop over element
243# include "vectorize.inc"
244 DO i = 1,nel
245 var_reg(i,k) = third*(dnl(pos1(i)+k-1)
246 . + dnl(pos2(i)+k-1)
247 . + dnl(pos3(i)+k-1))
248 ENDDO
249 ENDDO
250c
251 CALL c3evec3(1 ,nel ,surf,
252 . x1 ,x2 ,x3 ,y1 ,y2 ,
253 . y3 ,z1 ,z2 ,z3 ,
254 . e1x ,e2x ,e3x ,e1y ,e2y ,e3y ,
255 . e1z ,e2z ,e3z ,
256 . x31 ,y31 ,z31 ,x2l ,x3l ,y3l )
257C
258 ! Filling internal variable data of the non-local material
259 CALL cnloc_matini(elbuf_tab(ng),nel ,ipm ,
260 . bufmat ,time ,var_reg ,
261 . failure )
262C
263 !-----------------------------------------------------------------------
264 ! Computation of the element volume and the BtB matrix product
265 !-----------------------------------------------------------------------
266 ! Loop over elements
267# include "vectorize.inc"
268 DO i=1,nel
269c
270 ! Computation of the shape functions derivatives
271 px1(i) = -half*y3l(i)
272 py1(i) = half*(x3l(i)-x2l(i))
273 py2(i) = -half*x3l(i)
274c
275 ! Computation of the product BtxB
276 btb11(i) = px1(i)**2 + py1(i)**2
277 btb12(i) = -px1(i)**2 + py1(i)*py2(i)
278 btb13(i) = -py1(i)*(py1(i)+py2(i))
279 btb22(i) = px1(i)**2 + py2(i)**2
280 btb23(i) = -py2(i)*(py1(i)+py2(i))
281 btb33(i) = (py1(i)+py2(i))**2
282c
283 ! Computation of the element volume
284 vols(i) = area(numelc+nft+i)*thck(i)
285c
286 ! To check if element is not broken
287 offg(i) = elbuf_tab(ng)%GBUF%OFF(i)
288 ENDDO
289C
290 !-----------------------------------------------------------------------
291 ! Pre-treatment non-local regularization in the thickness
292 !-----------------------------------------------------------------------
293 ! Only if NDOF > 1
294 IF ((ndof > 1).AND.(len>zero)) THEN
295c
296 ! Pointing the non-local values in the thickness of the corresponding element
297 bufnl => elbuf_tab(ng)%NLOC(1,1)
298c
299 ! Pointing the non-local values in the thickness of the corresponding element
300 massth => bufnl%MASSTH(1:nel,1:ndnod)
301 unlth => bufnl%UNLTH(1:nel ,1:ndnod)
302 vnlth => bufnl%VNLTH(1:nel ,1:ndnod)
303 fnlth => bufnl%FNLTH(1:nel ,1:ndnod)
304c
305 DO k = 1,ndnod
306 DO i = 1,nel
307 ! Prediction of the velocities
308 vpred(i,k) = vnlth(i,k) - (fnlth(i,k)/massth(i,k))*(dt_nl/two)
309 ENDDO
310 ENDDO
311 DO k = 1,ndnod
312 DO i = 1,nel
313 ! Resetting non-local forces
314 fnlth(i,k) = zero
315 ENDDO
316 ENDDO
317c
318 ! Computation of non-local forces in the shell thickness
319 DO k = 1, ndof
320c
321 ! Computation of shape functions value
322 IF ((ndof==2).AND.(k==2)) THEN
323 nth1 = (z01(k,ndof) - zn1(k,ndof))/(zn1(k-1,ndof) - zn1(k,ndof))
324 nth2 = (z01(k,ndof) - zn1(k-1,ndof)) /(zn1(k,ndof) - zn1(k-1,ndof))
325 ELSE
326 nth1 = (z01(k,ndof) - zn1(k+1,ndof))/(zn1(k,ndof) - zn1(k+1,ndof))
327 nth2 = (z01(k,ndof) - zn1(k,ndof)) /(zn1(k+1,ndof) - zn1(k,ndof))
328 ENDIF
329c
330 ! Loop over elements
331 DO i = 1,nel
332 ! Computation of B-matrix values
333 IF ((ndof==2).AND.(k==2)) THEN
334 bth1 = (one/(zn1(k-1,ndof) - zn1(k,ndof)))*(one/thck(i))
335 bth2 = (one/(zn1(k,ndof) - zn1(k-1,ndof)))*(one/thck(i))
336 ELSE
337 bth1 = (one/(zn1(k,ndof) - zn1(k+1,ndof)))*(one/thck(i))
338 bth2 = (one/(zn1(k+1,ndof) - zn1(k,ndof)))*(one/thck(i))
339 ENDIF
340c
341 ! Computation of the non-local K matrix
342 k1 = (len**2)*(bth1**2) + nth1**2
343 k12 = (len**2)*(bth1*bth2)+ (nth1*nth2)
344 k2 = (len**2)*(bth2**2) + nth2**2
345c
346 ! Computation of the non-local forces
347 IF ((ndof==2).AND.(k==2)) THEN
348 fnlth(i,k-1) = fnlth(i,k-1) + (k1*unlth(i,k-1) + k12*unlth(i,k)
349 . + damp*((nth1**2)*vpred(i,k-1)
350 . + (nth1*nth2)*vpred(i,k))
351 . - (nth1*var_reg(i,k)))*vols(i)*wf1(k,ndof)
352 fnlth(i,k) = fnlth(i,k) + (k12*unlth(i,k-1) + k2*unlth(i,k)
353 . + damp*(nth1*nth2*vpred(i,k-1)
354 . + (nth2**2)*vpred(i,k))
355 . - nth2*var_reg(i,k))*vols(i)*wf1(k,ndof)
356 ELSE
357 fnlth(i,k) = fnlth(i,k) + (k1*unlth(i,k) + k12*unlth(i,k+1)
358 . + damp*((nth1**2)*vpred(i,k)
359 . + (nth1*nth2)*vpred(i,k+1))
360 . - (nth1*var_reg(i,k)))*vols(i)*wf1(k,ndof)
361 fnlth(i,k+1) = fnlth(i,k+1) + (k12*unlth(i,k) + k2*unlth(i,k+1)
362 . + damp*(nth1*nth2*vpred(i,k)
363 . + (nth2**2)*vpred(i,k+1))
364 . - nth2*var_reg(i,k))*vols(i)*wf1(k,ndof)
365 ENDIF
366 ENDDO
367 ENDDO
368c
369 DO k = 1,ndnod
370 DO i = 1,nel
371 ! Updating the non-local in-thickness velocities
372 vnlth(i,k) = vnlth(i,k) - (fnlth(i,k)/massth(i,k))*dt_nl
373 ENDDO
374 ENDDO
375c
376 DO k = 1,ndnod
377 DO i = 1,nel
378 ! Computing the non-local in-thickness cumulated values
379 unlth(i,k) = unlth(i,k) + vnlth(i,k)*dt_nl
380 ENDDO
381 ENDDO
382c
383 ! Transfer at the integration point
384 DO k = 1, ndof
385 !Computation of shape functions value
386 IF ((ndof==2).AND.(k==2)) THEN
387 nth1 = (z01(k,ndof) - zn1(k,ndof))/(zn1(k-1,ndof) - zn1(k,ndof))
388 nth2 = (z01(k,ndof) - zn1(k-1,ndof)) /(zn1(k,ndof) - zn1(k-1,ndof))
389 ELSE
390 nth1 = (z01(k,ndof) - zn1(k+1,ndof))/(zn1(k,ndof) - zn1(k+1,ndof))
391 nth2 = (z01(k,ndof) - zn1(k,ndof)) /(zn1(k+1,ndof) - zn1(k,ndof))
392 ENDIF
393 ! Loop over elements
394 DO i = 1,nel
395 !Integration points non-local variables
396 IF ((ndof==2).AND.(k==2)) THEN
397 var_reg(i,k) = nth1*unlth(i,k-1) + nth2*unlth(i,k)
398 ELSE
399 var_reg(i,k) = nth1*unlth(i,k) + nth2*unlth(i,k+1)
400 ENDIF
401 ENDDO
402 ENDDO
403 ENDIF
404c
405 !-----------------------------------------------------------------------
406 ! Computation of the elementary non-local forces
407 !-----------------------------------------------------------------------
408 ! Loop over integration points in the thickness
409 DO k = 1,ndof
410c
411 ! Computation of non-local forces
412# include "vectorize.inc"
413 DO i = 1,nel
414c
415 ! If the element is not broken
416 IF (offg(i) > zero) THEN
417 ! Computing the product BtB*UNL
418 b1 = ((len**2)/vols(i))*wf1(k,ndof)*(btb11(i)*unl(pos1(i)+k-1) + btb12(i)*unl(pos2(i)+k-1)
419 . + btb13(i)*unl(pos3(i)+k-1))
420 b2 = ((len**2)/vols(i))*wf1(k,ndof)*(btb12(i)*unl(pos1(i)+k-1) + btb22(i)*unl(pos2(i)+k-1)
421 . + btb23(i)*unl(pos3(i)+k-1))
422 b3 = ((len**2)/vols(i))*wf1(k,ndof)*(btb13(i)*unl(pos1(i)+k-1) + btb23(i)*unl(pos2(i)+k-1)
423 . + btb33(i)*unl(pos3(i)+k-1))
424 ! Computing the product NtN*UNL
425 ntn_unl = ((unl(pos1(i)+k-1) + unl(pos2(i)+k-1) + unl(pos3(i)+k-1))*third*third)*vols(i)*wf1(k,ndof)
426 ! Computing the product DAMP*NtN*VNL!
427 ntn_vnl = ((vnl(pos1(i)+k-1) + vnl(pos2(i)+k-1) + vnl(pos3(i)+k-1))*third*third)*damp*vols(i)*wf1(k,ndof)
428 ! Introducing the internal variable to be regularized
429 ntvar = var_reg(i,k)*third*vols(i)*wf1(k,ndof)
430 !Computing the elementary non-local forces
431 fnl(pos1(i)+k-1) = fnl(pos1(i)+k-1) - (ntn_unl + ntn_vnl - ntvar + b1)
432 fnl(pos2(i)+k-1) = fnl(pos2(i)+k-1) - (ntn_unl + ntn_vnl - ntvar + b2)
433 fnl(pos3(i)+k-1) = fnl(pos3(i)+k-1) - (ntn_unl + ntn_vnl - ntvar + b3)
434 ELSE
435 ! Non-local absorbing forces
436 fnl(pos1(i)+k-1) = fnl(pos1(i)+k-1) -
437 . wf1(k,ndof)*dens*sspnl*vnl(pos1(i)+k-1)*sqrt((four/sqrt(three))*(le_max**2))*thck(i)
438 fnl(pos2(i)+k-1) = fnl(pos2(i)+k-1) -
439 . wf1(k,ndof)*dens*sspnl*vnl(pos2(i)+k-1)*sqrt((four/sqrt(three))*(le_max**2))*thck(i)
440 fnl(pos3(i)+k-1) = fnl(pos3(i)+k-1) -
441 . wf1(k,ndof)*dens*sspnl*vnl(pos3(i)+k-1)*sqrt((four/sqrt(three))*(le_max**2))*thck(i)
442 ENDIF
443 ENDDO
444 ENDDO
445c -------------------
446 IF (ALLOCATED(var_reg)) DEALLOCATE(var_reg)
447 IF (ALLOCATED(vpred)) DEALLOCATE(vpred)
subroutine cnloc_matini(elbuf_str, nel, ipm, bufmat, time, varnl, failure)
#define my_real
Definition cppsort.cpp:32
subroutine area(d1, x, x2, y, y2, eint, stif0)
#define min(a, b)
Definition macros.h:20
subroutine c3evec3(jft, jlt, area, x1, x2, x3, y1, y2, y3, z1, z2, z3, e1x, e2x, e3x, e1y, e2y, e3y, e1z, e2z, e3z, x31, y31, z31, x2l, x3l, y3l)
Definition c3evec3.F:39