Rocstar  1.0
Rocstar multiphysics simulation application
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
v3d4_ale.f90
Go to the documentation of this file.
1 !*********************************************************************
2 !* Illinois Open Source License *
3 !* *
4 !* University of Illinois/NCSA *
5 !* Open Source License *
6 !* *
7 !* Copyright@2008, University of Illinois. All rights reserved. *
8 !* *
9 !* Developed by: *
10 !* *
11 !* Center for Simulation of Advanced Rockets *
12 !* *
13 !* University of Illinois *
14 !* *
15 !* www.csar.uiuc.edu *
16 !* *
17 !* Permission is hereby granted, free of charge, to any person *
18 !* obtaining a copy of this software and associated documentation *
19 !* files (the "Software"), to deal with the Software without *
20 !* restriction, including without limitation the rights to use, *
21 !* copy, modify, merge, publish, distribute, sublicense, and/or *
22 !* sell copies of the Software, and to permit persons to whom the *
23 !* Software is furnished to do so, subject to the following *
24 !* conditions: *
25 !* *
26 !* *
27 !* @ Redistributions of source code must retain the above copyright *
28 !* notice, this list of conditions and the following disclaimers. *
29 !* *
30 !* @ Redistributions in binary form must reproduce the above *
31 !* copyright notice, this list of conditions and the following *
32 !* disclaimers in the documentation and/or other materials *
33 !* provided with the distribution. *
34 !* *
35 !* @ Neither the names of the Center for Simulation of Advanced *
36 !* Rockets, the University of Illinois, nor the names of its *
37 !* contributors may be used to endorse or promote products derived *
38 !* from this Software without specific prior written permission. *
39 !* *
40 !* THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, *
41 !* EXPRESS OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES *
42 !* OF MERCHANTABILITY, FITNESS FOR A PARTICULAR PURPOSE AND *
43 !* NONINFRINGEMENT. IN NO EVENT SHALL THE CONTRIBUTORS OR *
44 !* COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER *
45 !* LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, *
46 !* ARISING FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE *
47 !* USE OR OTHER DEALINGS WITH THE SOFTWARE. *
48 !*********************************************************************
49 !* Please acknowledge The University of Illinois Center for *
50 !* Simulation of Advanced Rockets in works and publications *
51 !* resulting from this software or its derivatives. *
52 !*********************************************************************
53 SUBROUTINE v3d4_ale(v_bar,a_bar,d,vhalf,Rnet, &
54  e,xnu,rho,numnp,numat_vol, &
55  numcstet,matcstet,lmcstet,meshcoor, &
56  nstart,nend)
57 
58  IMPLICIT NONE
59 
60  INTEGER :: numnp ! number of nodes
61  INTEGER :: numat_vol ! number of materials
62  INTEGER :: numcstet ! number of tetrahedrals
63 ! nodal mesh velocities
64  REAL*8, DIMENSION(1:3*numnp) :: v_bar
65 ! nodal mesh accelerations
66  REAL*8, DIMENSION(1:3*numnp) :: a_bar
67 ! nodal displacements
68  REAL*8, DIMENSION(1:3*numnp) :: d
69 ! nodal velocities
70  REAL*8, DIMENSION(1:3*numnp) :: vhalf
71 ! Contibution to the Net Force Vector
72  REAL*8, DIMENSION(1:3*numnp) :: rnet
73 ! young's moduli
74  REAL*8, DIMENSION(1:numat_vol) :: e
75 ! Poisson's ratios
76  REAL*8, DIMENSION(1:numat_vol) :: xnu
77 ! densities
78  REAL*8, DIMENSION(1:numat_vol) :: rho
79 ! mesh coordinates
80  REAL*8, DIMENSION(1:3,1:numnp) :: meshcoor
81 ! Tet connectivity table
82  INTEGER, DIMENSION(1:4,1:numcstet) :: lmcstet
83 
84  INTEGER matcstet(numcstet) ! list of all materials
85  INTEGER ielem , nx, ny, nz
86 
87 ! CODE FOR 3-D ALE, ELASTODYNAMIC, FIRST-ORDER SOLID ELEMENT
88 
89 ! integer ndof, nnode, ndim, nmat, nsurf
90  INTEGER ndof, nnode, ndim, nmat, nsurf
91 
92 !-- for tetrahedron element !
93  REAL*8 meshpos(3,4), meshvel(3,4), meshacc(3,4), mat(3)
94 
95 ! INPUT:
96 
97 ! ndof -- # of element degrees of freedom
98 ! nnode -- # of nodes on element
99 ! ndim -- # of spatial dimensions
100 ! nmat -- # of material parameters (rho, E, nu)
101 ! nsurf -- # of surfaces on element
102 ! meshpos -- coordinates of position w.r.t. to global basis of nodes
103 ! meshvel -- components of mesh velocity w.r.t to global basis at nodes
104 ! meshacc -- '' '' '' accelern. '' '' '' '' '' ''
105 ! mat -- array containing material parameters (e.g. rho, E, nu)
106 
107  INTEGER nintk_v, knode, &
108  knodr, knodc,i,j,rindx,cindx,k,mm,ll, &
109  iface, igauss, grindx,grindy,grindz
110 !
111 !--- new variables for tetrahedron
112  parameter( nintk_v = 1) ! 4-noded tetrahedron
113 
114  DOUBLE PRECISION shfn(8,8), shdx(3,8,8), dv(8), dvsum, &
115  mvgss(3), magss(3), mcgss(3,3), &
116  wsp1,wsparr(3),trmcgss,wsp2,wsp3,wsp4, &
117  third,wsp5,vxi(8),veta(8),vzeta(8),wvar, &
118  sxi(4,6),seta(4,6),szeta(4,6),shdxi(3,8,8), &
119  wsparr1(3),shdx_av(3,8)
120  DOUBLE PRECISION disp(12), velo(12), rimsi(12), dimsi(12)
121  DOUBLE PRECISION multiplier , ajacin(3,3), surfnormal(3,4)
122  DOUBLE PRECISION wt_cst, da, sum,factor
123  DOUBLE PRECISION sto_wsp1(4),sto_wsp2(4),sto_wsp3(4),sto_wsp4(4)
124  INTEGER surf_elem(3,4), k2
125 
126  INTEGER :: nstart, nend
127 
128  ndim = 3
129  nnode = 4
130  ndof = ndim * nnode
131  nsurf = 4
132  multiplier = 1.d0/6.d0
133 
134 
135 ! initialize integration pt.coords
136 ! one point integration rule for tetrahedron
137 ! reference: The Finite Element Method by Thomas J.R. Hughes
138 ! Chap.3 pp 170-174
139 
140  vxi(1) = 0.25d0
141  veta(1) = 0.25d0
142  vzeta(1) = 0.25d0
143 
144 !-- one point on each surface
145 ! face 1 : node 1 2 3
146 ! face 2 : node 1 2 4
147 ! face 3 : node 2 3 4
148 ! face 4 : node 1 3 4
149 
150  wvar = 1.d0/dsqrt(3.d0)
151  third = 1.d0/3.d0
152 
153  sxi(1,1) = third
154  seta(1,1) = 0.d0
155  szeta(1,1) = third
156 
157  sxi(1,2) = third
158  seta(1,2) = third
159  szeta(1,2) = third
160 
161  sxi(1,3) = third
162  seta(1,3) = third
163  szeta(1,3) = 0.d0
164 
165  sxi(1,4) = 0.d0
166  seta(1,4) = third
167  szeta(1,4) = third
168 
169 !-- surface normal vectors
170  surf_elem(1,1)=1
171  surf_elem(2,1)=3
172  surf_elem(3,1)=2
173 
174  surf_elem(1,2)=1
175  surf_elem(2,2)=2
176  surf_elem(3,2)=4
177 
178  surf_elem(1,3)=3
179  surf_elem(2,3)=4
180  surf_elem(3,3)=2
181 
182  surf_elem(1,4)=1
183  surf_elem(2,4)=4
184  surf_elem(3,4)=3
185 
186 
187 ! loop for elements
188  igauss = 1
189  DO ielem = nstart, nend
190  j = matcstet(ielem)
191  mat(1)= rho(j)
192  mat(2)= e(j)
193  mat(3)= xnu(j)
194 
195  DO i=1, nnode
196  nx = lmcstet(i,ielem)*3-2
197  ny = lmcstet(i,ielem)*3-1
198  nz = lmcstet(i,ielem)*3
199  meshpos(1,i) = meshcoor(1,lmcstet(i,ielem))
200  meshpos(2,i) = meshcoor(2,lmcstet(i,ielem))
201  meshpos(3,i) = meshcoor(3,lmcstet(i,ielem))
202  meshvel(1,i) = v_bar(nx)
203  meshvel(2,i) = v_bar(ny)
204  meshvel(3,i) = v_bar(nz)
205  meshacc(1,i) = a_bar(nx)
206  meshacc(2,i) = a_bar(ny)
207  meshacc(3,i) = a_bar(nz)
208  disp(i*3-2) = d(nx)
209  disp(i*3-1) = d(ny)
210  disp(i*3 ) = d(nz)
211  velo(i*3-2) = vhalf(nx)
212  velo(i*3-1) = vhalf(ny)
213  velo(i*3 ) = vhalf(nz)
214  END DO
215 
216 
217  CALL shcalc(shfn, shdx, dv, dvsum, shdx_av, meshpos, &
218  ndim,nnode,nintk_v,vxi,veta,vzeta,shdxi,ajacin)
219 
220 
221 ! fill in element matrices for contributions from
222 ! volume integration points.
223 
224 
225 ! evaluate mesh-related quantitites at vol. int. pt.
226 
227  mvgss(1) = 0.25d0*(meshvel(1,1)+meshvel(1,2)+ &
228  meshvel(1,3)+meshvel(1,4))
229  mvgss(2) = 0.25d0*(meshvel(2,1)+meshvel(2,2)+ &
230  meshvel(2,3)+meshvel(2,4))
231  mvgss(3) = 0.25d0*(meshvel(3,1)+meshvel(3,2)+ &
232  meshvel(3,3)+meshvel(3,4))
233  magss(1) = 0.25d0*(meshacc(1,1)+meshacc(1,2)+ &
234  meshacc(1,3)+meshacc(1,4))
235  magss(2) = 0.25d0*(meshacc(2,1)+meshacc(2,2)+ &
236  meshacc(2,3)+meshacc(2,4))
237  magss(3) = 0.25d0*(meshacc(3,1)+meshacc(3,2)+ &
238  meshacc(3,3)+meshacc(3,4))
239 
240  mcgss(1,1) = 0.d0
241  mcgss(1,2) = 0.d0
242  mcgss(1,3) = 0.d0
243  mcgss(2,1) = 0.d0
244  mcgss(2,2) = 0.d0
245  mcgss(2,3) = 0.d0
246  mcgss(3,1) = 0.d0
247  mcgss(3,2) = 0.d0
248  mcgss(3,3) = 0.d0
249  DO knode = 1,nnode
250  mcgss(1,1) = mcgss(1,1) + meshvel(1,knode)* &
251  shdx(1,knode,igauss)
252  mcgss(2,1) = mcgss(2,1) + meshvel(2,knode)* &
253  shdx(1,knode,igauss)
254  mcgss(3,1) = mcgss(3,1) + meshvel(3,knode)* &
255  shdx(1,knode,igauss)
256  mcgss(1,2) = mcgss(1,2) + meshvel(1,knode)* &
257  shdx(2,knode,igauss)
258  mcgss(2,2) = mcgss(2,2) + meshvel(2,knode)* &
259  shdx(2,knode,igauss)
260  mcgss(3,2) = mcgss(3,2) + meshvel(3,knode)* &
261  shdx(2,knode,igauss)
262  mcgss(1,3) = mcgss(1,3) + meshvel(1,knode)* &
263  shdx(3,knode,igauss)
264  mcgss(2,3) = mcgss(2,3) + meshvel(2,knode)* &
265  shdx(3,knode,igauss)
266  mcgss(3,3) = mcgss(3,3) + meshvel(3,knode)* &
267  shdx(3,knode,igauss)
268  END DO
269 
270  trmcgss = mcgss(1,1) + mcgss(2,2) + mcgss(3,3)
271 
272  wsparr(1) = mcgss(1,1)*mvgss(1) + &
273  mcgss(1,2)*mvgss(2) + &
274  mcgss(1,3)*mvgss(3)
275  wsparr(2) = mcgss(2,1)*mvgss(1) + &
276  mcgss(2,2)*mvgss(2) + &
277  mcgss(2,3)*mvgss(3)
278  wsparr(3) = mcgss(3,1)*mvgss(1) + &
279  mcgss(3,2)*mvgss(2) + &
280  mcgss(3,3)*mvgss(3)
281 
282 
283  DO knodc = 1,nnode
284  sto_wsp1(knodc) = &
285  mvgss(1)*shdx(1,knodc,igauss) + &
286  mvgss(2)*shdx(2,knodc,igauss) + &
287  mvgss(3)*shdx(3,knodc,igauss)
288  sto_wsp2(knodc) = &
289  wsparr(1)*shdx(1,knodc,igauss) + &
290  wsparr(2)*shdx(2,knodc,igauss) + &
291  wsparr(3)*shdx(3,knodc,igauss)
292  sto_wsp3(knodc) = &
293  mvgss(1)*shdx(1,knodc,igauss) + &
294  mvgss(2)*shdx(2,knodc,igauss) + &
295  mvgss(3)*shdx(3,knodc,igauss)
296  sto_wsp4(knodc) = &
297  magss(1)*shdx(1,knodc,igauss) + &
298  magss(2)*shdx(2,knodc,igauss) + &
299  magss(3)*shdx(3,knodc,igauss)
300  END DO
301 
302  factor=mat(1)*dv(1)*multiplier
303  DO knodr = 1,nnode
304  grindx = (lmcstet(knodr,ielem)-1)*3 + 1
305  grindy = (lmcstet(knodr,ielem)-1)*3 + 2
306  grindz = (lmcstet(knodr,ielem)-1)*3 + 3
307 
308 ! damping
309 ! ----------------------------------
310 
311  rnet(grindx) = rnet(grindx) + &
312  factor*shfn(knodr,igauss)*2.d0* &
313  ( sto_wsp1(1)*velo((1-1)*ndim+1) + &
314  sto_wsp1(2)*velo((2-1)*ndim+1) + &
315  sto_wsp1(3)*velo((3-1)*ndim+1) + &
316  sto_wsp1(4)*velo((4-1)*ndim+1) )
317 
318  rnet(grindy) = rnet(grindy) + &
319  factor*shfn(knodr,igauss)*2.d0* &
320  ( sto_wsp1(1)*velo((1-1)*ndim+2) + &
321  sto_wsp1(2)*velo((2-1)*ndim+2) + &
322  sto_wsp1(3)*velo((3-1)*ndim+2) + &
323  sto_wsp1(4)*velo((4-1)*ndim+2) )
324 
325  rnet(grindz) = rnet(grindz) + &
326  factor*shfn(knodr,igauss)*2.d0* &
327  ( sto_wsp1(1)*velo((1-1)*ndim+3) + &
328  sto_wsp1(2)*velo((2-1)*ndim+3) + &
329  sto_wsp1(3)*velo((3-1)*ndim+3) + &
330  sto_wsp1(4)*velo((4-1)*ndim+3) )
331 
332 ! stiffness (vol. integration terms)
333 ! ----------------------------------
334 
335  rnet(grindx) = rnet(grindx) &
336  - factor*shfn(knodr,igauss) * ( &
337  -trmcgss*sto_wsp1(1)*disp((1-1)*ndim+1) &
338  -trmcgss*sto_wsp1(2)*disp((2-1)*ndim+1) &
339  -trmcgss*sto_wsp1(3)*disp((3-1)*ndim+1) &
340  -trmcgss*sto_wsp1(4)*disp((4-1)*ndim+1) &
341  + sto_wsp2(1)*disp((1-1)*ndim+1) &
342  + sto_wsp2(2)*disp((2-1)*ndim+1) &
343  + sto_wsp2(3)*disp((3-1)*ndim+1) &
344  + sto_wsp2(4)*disp((4-1)*ndim+1) &
345  - sto_wsp4(1)*disp((1-1)*ndim+1) &
346  - sto_wsp4(2)*disp((2-1)*ndim+1) &
347  - sto_wsp4(3)*disp((3-1)*ndim+1) &
348  - sto_wsp4(4)*disp((4-1)*ndim+1) ) &
349  - factor*sto_wsp3(knodr) * ( &
350  - sto_wsp1(1)*disp((1-1)*ndim+1) &
351  - sto_wsp1(2)*disp((2-1)*ndim+1) &
352  - sto_wsp1(3)*disp((3-1)*ndim+1) &
353  - sto_wsp1(4)*disp((4-1)*ndim+1) )
354 
355 
356  rnet(grindy) = rnet(grindy) &
357  - factor*shfn(knodr,igauss) * ( &
358  -trmcgss*sto_wsp1(1)*disp((1-1)*ndim+2) &
359  -trmcgss*sto_wsp1(2)*disp((2-1)*ndim+2) &
360  -trmcgss*sto_wsp1(3)*disp((3-1)*ndim+2) &
361  -trmcgss*sto_wsp1(4)*disp((4-1)*ndim+2) &
362  + sto_wsp2(1)*disp((1-1)*ndim+2) &
363  + sto_wsp2(2)*disp((2-1)*ndim+2) &
364  + sto_wsp2(3)*disp((3-1)*ndim+2) &
365  + sto_wsp2(4)*disp((4-1)*ndim+2) &
366  - sto_wsp4(1)*disp((1-1)*ndim+2) &
367  - sto_wsp4(2)*disp((2-1)*ndim+2) &
368  - sto_wsp4(3)*disp((3-1)*ndim+2) &
369  - sto_wsp4(4)*disp((4-1)*ndim+2) ) &
370  - factor*sto_wsp3(knodr) * ( &
371  - sto_wsp1(1)*disp((1-1)*ndim+2) &
372  - sto_wsp1(2)*disp((2-1)*ndim+2) &
373  - sto_wsp1(3)*disp((3-1)*ndim+2) &
374  - sto_wsp1(4)*disp((4-1)*ndim+2) )
375 
376  rnet(grindz) = rnet(grindz) &
377  - factor*shfn(knodr,igauss) * ( &
378  -trmcgss*sto_wsp1(1)*disp((1-1)*ndim+3) &
379  -trmcgss*sto_wsp1(2)*disp((2-1)*ndim+3) &
380  -trmcgss*sto_wsp1(3)*disp((3-1)*ndim+3) &
381  -trmcgss*sto_wsp1(4)*disp((4-1)*ndim+3) &
382  + sto_wsp2(1)*disp((1-1)*ndim+3) &
383  + sto_wsp2(2)*disp((2-1)*ndim+3) &
384  + sto_wsp2(3)*disp((3-1)*ndim+3) &
385  + sto_wsp2(4)*disp((4-1)*ndim+3) &
386  - sto_wsp4(1)*disp((1-1)*ndim+3) &
387  - sto_wsp4(2)*disp((2-1)*ndim+3) &
388  - sto_wsp4(3)*disp((3-1)*ndim+3) &
389  - sto_wsp4(4)*disp((4-1)*ndim+3) ) &
390  - factor*sto_wsp3(knodr) * ( &
391  - sto_wsp1(1)*disp((1-1)*ndim+3) &
392  - sto_wsp1(2)*disp((2-1)*ndim+3) &
393  - sto_wsp1(3)*disp((3-1)*ndim+3) &
394  - sto_wsp1(4)*disp((4-1)*ndim+3) )
395 
396 ! --------------------------------------------
397 ! all stiffness terms above are from physical
398 ! acceleration. The isotropic linear
399 ! elastic stiffness is calculated in CST
400 ! --------------------------------------------
401 
402  END DO
403 
404 
405 ! contribution to stiffness matrix from surface terms
406 ! ---------------------------------------------------
407 
408  DO iface = 1,nsurf
409 
410 !-- shdx is constant and already calculated !!
411 
412  shfn(1,1)=szeta(1,iface)
413  shfn(2,1)=sxi(1,iface)
414  shfn(3,1)=1.d0-sxi(1,iface)-seta(1,iface)-szeta(1,iface)
415  shfn(4,1)=seta(1,iface)
416 
417  DO j = 1,ndim
418  mvgss(j) = 0.d0
419  DO knode = 1,nnode
420  mvgss(j) = mvgss(j) + meshvel(j,knode)* &
421  shfn(knode,igauss)
422  END DO
423  END DO
424 
425 !-- unit out-normal vector on each surface
426 
427  DO j = 1, 3
428  wsparr(j) =meshpos(j,surf_elem(2,iface)) - &
429  meshpos(j,surf_elem(1,iface))
430  wsparr1(j)=meshpos(j,surf_elem(3,iface)) - &
431  meshpos(j,surf_elem(1,iface))
432  END DO
433 
434 ! do cross product: wsparr = nk * da
435 
436  wsp1 = wsparr(2)*wsparr1(3) - wsparr(3)*wsparr1(2)
437  wsp2 = wsparr(3)*wsparr1(1) - wsparr(1)*wsparr1(3)
438  wsp3 = wsparr(1)*wsparr1(2) - wsparr(2)*wsparr1(1)
439  wsparr(1) = wsp1 * 0.5d0
440  wsparr(2) = wsp2 * 0.5d0
441  wsparr(3) = wsp3 * 0.5d0
442 
443 
444 ! evaluate terms required in stiffness
445 
446  wsp1 = mvgss(1)*wsparr(1) &
447  + mvgss(2)*wsparr(2) &
448  + mvgss(3)*wsparr(3)
449 
450 ! fill stiffness matrix
451 
452 
453  DO knodc = 1,nnode
454  sto_wsp2(knodc) = &
455  mvgss(1)*shdx(1,knodc,igauss) + &
456  mvgss(2)*shdx(2,knodc,igauss) + &
457  mvgss(3)*shdx(3,knodc,igauss)
458  END DO
459 
460 
461  DO knodr = 1, nnode
462  grindx = (lmcstet(knodr,ielem)-1)*3 + 1
463  grindy = (lmcstet(knodr,ielem)-1)*3 + 2
464  grindz = (lmcstet(knodr,ielem)-1)*3 + 3
465 
466  rnet(grindx) = rnet(grindx) &
467  - mat(1)*shfn(knodr,igauss)*wsp1*( &
468  + sto_wsp2(1)*disp((1-1)*ndim+1) &
469  + sto_wsp2(2)*disp((2-1)*ndim+1) &
470  + sto_wsp2(3)*disp((3-1)*ndim+1) &
471  + sto_wsp2(4)*disp((4-1)*ndim+1) )
472 
473  rnet(grindy) = rnet(grindy) &
474  - mat(1)*shfn(knodr,igauss)*wsp1*( &
475  + sto_wsp2(1)*disp((1-1)*ndim+2) &
476  + sto_wsp2(2)*disp((2-1)*ndim+2) &
477  + sto_wsp2(3)*disp((3-1)*ndim+2) &
478  + sto_wsp2(4)*disp((4-1)*ndim+2) )
479 
480  rnet(grindz) = rnet(grindz) &
481  - mat(1)*shfn(knodr,igauss)*wsp1*( &
482  + sto_wsp2(1)*disp((1-1)*ndim+3) &
483  + sto_wsp2(2)*disp((2-1)*ndim+3) &
484  + sto_wsp2(3)*disp((3-1)*ndim+3) &
485  + sto_wsp2(4)*disp((4-1)*ndim+3) )
486  END DO
487 
488  END DO ! for iface=1,nsurf
489 
490 1112 CONTINUE
491 
492  END DO ! for ielem=1,numcstet
493 
494 111 FORMAT(i5,5e15.7)
495 
496  RETURN
497 END SUBROUTINE v3d4_ale
498 
Tfloat sum() const
Return the sum of all the pixel values in an image.
Definition: CImg.h:13022
const NT & d
j indices k indices k
Definition: Indexing.h:6
blockLoc i
Definition: read.cpp:79
subroutine shcalc(shfn, shdx, dv, dvsum, shdx_av, meshpos, ndim, nnode, nintk, xi, eta, zeta, shdxi, wsp2)
Definition: shcalc.f90:53
j indices j
Definition: Indexing.h:6
subroutine v3d4_ale(v_bar, a_bar, d, vhalf, Rnet, E, xnu, rho, numnp, numat_vol, numcstet, matcstet, lmcstet, meshcoor, nstart, nend)
Definition: v3d4_ale.f90:53