V3FIT
All Classes Namespaces Files Functions Variables Enumerations Macros Pages
smpint.f
1  SUBROUTINE smpint( ND, NF, MINCLS, MAXCLS, FUNSUB,
2  & EPSABS, EPSREL, KEY, SBRGNS, WRKLEN, VRTWRK,
3  & RESTAR, VALUE, ERROR, FUNCLS, INFORM )
4 *
5 ****BEGIN PROLOGUE SMPINT
6 ****AUTHOR
7 *
8 * Alan Genz
9 * Department of Mathematics
10 * Washington State University
11 * Pullman, WA 99164-3113, USA
12 * Email: alangenz@wsu.edu
13 *
14 ****LAST MODIFICATION 2002-9
15 ****KEYWORDS automatic multidimensional integrator,
16 * n-dimensional simplex, general purpose, global adaptive
17 ****PURPOSE To calculate an approximation to a vector of integrals
18 *
19 * I = I (F ,F ,...,F ) DS
20 * S 1 2 NF
21 *
22 * where S is a collection ND-dimensional simplices,
23 * and F = F (X ,X ,...,X ), K = 1, 2, ..., NF,
24 * K K 1 2 ND
25 * and try to satisfy for each component I(K) of I
26 * ABS( I(K) - VALUE(K) ) < MAX( EPSABS, EPSREL*ABS(I(K)) )
27 *
28 ****DESCRIPTION Computation of integrals over simplical regions.
29 * SMPINT is a driver for the integration routine SMPSAD,
30 * which repeatedly subdivides the region of integration and
31 * estimates the integrals and the errors over the subregions
32 * with greatest estimated errors until the error request
33 * is met or MAXCLS function evaluations have been used.
34 *
35 * ON ENTRY
36 *
37 * ND Integer, number of variables. 1 < ND
38 * NF Integer, number of components of the integral.
39 * MINCLS Integer, minimum number of FUNSUB calls.
40 * MAXCLS Integer, maximum number of FUNSUB calls.
41 * RULCLS is number FUNSUB calls for each subregion (see WRKLEN),
42 * DIFCLS = 1 + 2*ND*( ND + 1 ).
43 * If RESTAR = 0, MAXCLS must be >= MAX(SBRGNS*RULCLS,MINCLS).
44 * If RESTAR = 1, MAXCLS must be >= MAX(4*RULCLS+DIFCLS,MINCLS).
45 * FUNSUB Externally declared subroutine for computing components of
46 * the integrand at the given evaluation point.
47 * It must have parameters (ND,X,NF,FUNVLS)
48 * Input parameters:
49 * ND Integer that gives the dimension of I
50 * X Real array of dimension ND that contains the
51 * evaluation point.
52 * NF Integer that gives the number of components of I.
53 * Output parameter:
54 * FUNVLS Real array of dimension NF that contains the
55 * components of the integrand.
56 * EPSABS Real.
57 * Requested absolute accuracy.
58 * EPSREL Real requested relative accuracy.
59 * KEY Integer, key to selected local integration rule.
60 * KEY = 3 gives the user a (default) degree 7 integration rule.
61 * KEY = 1 gives the user a degree 3 integration rule.
62 * KEY = 2 gives the user a degree 5 integration rule.
63 * KEY = 3 gives the user a degree 7 integration rule.
64 * KEY = 4 gives the user a degree 9 integration rule.
65 * WRKLEN Integer, length of the working array VRTWRK.
66 * WRKLEN should be >= WRKSBS*( ND*(ND+1) + 2*NF + 3 )
67 * + (ND+1)*(ND+2) + 7*NF, where
68 * WRKSBS = SBRGNS + 3*( MAXCLS/RULCLS - SBRGNS*(1-RESTAR) )/4.
69 * If
70 * KEY = 0, RULCLS = (ND+4)*(ND+3)*(ND+2)/6 + (ND+2)*(ND+1);
71 * KEY = 1, RULCLS = 2*ND+3;
72 * KEY = 2, RULCLS = (ND+3)*(ND+2)/2 + 2*(ND+1);
73 * KEY = 3, RULCLS = (ND+4)*(ND+3)*(ND+2)/6 + (ND+2)*(ND+1);
74 * KEY = 4, RULCLS = (ND+5)*(ND+4)*(ND+3)*(ND+2)/24
75 * + 5*(ND+2)*(ND+1)/2 .
76 * SBRGNS Integer, initial number of simplices.
77 * VRTWRK Real array of dimension WRKLEN.
78 * Work should contain the simplex vertices for SBRGNS
79 * simplices; the coordinates of vertex J for simplex K
80 * must be in VRTWRK(I+J*ND+(K-1)*ND*(ND+1)),
81 * for I = 1, ..., ND; J = 0, ..., ND; K = 1, ..., SBRGNS.
82 * The rest of VRTWRK is used as working storage; see below.
83 * RESTAR Integer.
84 * If RESTAR = 0, this is the first attempt to compute
85 * the integrals over SBRGNS simplices.
86 * If RESTAR = 1, then we restart a previous attempt.
87 * In this case the only parameters for SMPINT that may
88 * be changed (with respect to the previous call of SMPINT)
89 * are MINCLS, MAXCLS, EPSABS, EPSREL, KEY and RESTAR.
90 * ON RETURN
91 *
92 * SBRGNS Integer.
93 * SBRGNS contains the current number of simplices. They
94 * were obtained by subdividing the input simplicies.
95 * VRTWRK Real array of dimension WRKLEN.
96 * Used as working storage.
97 * Let MAXSUB = (WRKLEN-(ND+1)*(ND+2)-7*NF)/(ND*(ND+1)+2*NF+3).
98 * VRTWRK(1), ..., VRTWRK(ND*(ND+1)), ...,
99 * VRTWRK(ND*(ND+1)*SBRGNS) contain subregion vertices.
100 * VRTWRK(ND*(ND+1)*MAXSUB+1), ...,
101 * VRTWRK(ND*(ND+1)*MAXSUB+NF*SBRGNS) contain
102 * estimated components of the integrals over the subregions.
103 * VRTWRK(ND*(ND+1)*MAXSUB+NF*MAXSUB+1), ...,
104 * VRTWRK(ND*(ND+1)*MAXSUB+NF*(MAXSUB+SBRGNS))
105 * contain estimated errors for the subregions.
106 * VRTWRK(ND*(ND+1)*MAXSUB+2*NF*MAXSUB+1), ...,
107 * VRTWRK(ND*(ND+1)*MAXSUB+2*NF*MAXSUB+SBRGNS))
108 * contain volumes of the subregions (scaled by ND!).
109 * VRTWRK(ND*(ND+1)*MAXSUB+(2*NF+1)*MAXSUB+1), ...,
110 * VRTWRK(ND*(ND+1)*MAXSUB+(2*NF+1)*MAXSUB+SBRGNS)
111 * contain greatest errors in each subregion.
112 * VRTWRK(ND*(ND+1)*MAXSUB+(2*NF+2)*MAXSUB+1), ...,
113 * VRTWRK(ND*(ND+1)*MAXSUB+(2*NF+2)*MAXSUB+SBRGNS)
114 * contain pointers for the subregion heap.
115 * The rest of VRTWRK is used as temporary storage in SMPSAD.
116 * VALUE Real array of dimension NF of integral approximations.
117 * ERROR Real array of dimension NF, of absolute accuracy estimates.
118 * FUNCLS Integer, number of FUNSUB calls used by SMPINT.
119 * INFORM Integer.
120 * INFORM = 0 for normal exit, when ERROR(K) <= EPSABS or
121 * ERROR(K) <= ABS(VALUE(K))*EPSREL with MAXCLS or less
122 * function evaluations for all values of K, 1 <= K <= NF.
123 * INFORM = 1 if MAXCLS was too small for SMPINT to obtain
124 * the required accuracy. In this case SMPINT returns
125 * values VALUE with estimated absolute accuracies ERROR.
126 * INFORM = 2 if KEY < 0 or KEY > 4,
127 * INFORM = 3 if ND < 2,
128 * INFORM = 4 if NF < 1,
129 * INFORM = 5 if EPSABS < 0 and EPSREL < 0,
130 * INFORM = 6 if WRKLEN is too small,
131 * INFORM = 7 if RESTAR < 0 or RESTAR > 1,
132 * INFORM = 8 if SBRGNS <= 0.
133 *
134 ****ROUTINES CALLED SMPCHC,SMPSAD
135 ****END PROLOGUE SMPINT
136 *
137 * Global variables.
138 *
139  EXTERNAL funsub
140  INTEGER ND, NF, MINCLS, MAXCLS, SBRGNS
141  INTEGER KEY, WRKLEN, RESTAR, FUNCLS, INFORM
142  DOUBLE PRECISION EPSABS, EPSREL
143  DOUBLE PRECISION VALUE(NF), ERROR(NF), VRTWRK(WRKLEN)
144 *
145 * Local variables.
146 *
147 * MAXSUB Integer, maximum allowed number of subdivisions
148 * for the given values of KEY, ND and NF.
149 *
150  INTEGER MAXSUB, RULCLS, I1, I2, I3, I4, I5, I6
151 *
152 ****FIRST PROCESSING STATEMENT SMPINT
153 *
154 * Compute MAXSUB and RULCLS, and check the input parameters.
155 *
156 *
157  CALL smpchc( nd, nf, mincls, maxcls, epsabs, epsrel, sbrgns,
158  & key, wrklen, restar, rulcls, maxsub, inform )
159  IF ( inform .EQ. 0 ) THEN
160 *
161 * Split up the work space and call SMPSAD.
162 *
163  i1 = 1 + maxsub*nd*(nd+1)
164  i2 = i1 + maxsub*nf
165  i3 = i2 + maxsub*nf
166  i4 = i3 + maxsub
167  i5 = i4 + maxsub
168  i6 = i5 + maxsub
169  CALL smpsad( nd, nf, funsub, mincls, maxcls, epsabs, epsrel,
170  & restar, key, rulcls, maxsub, sbrgns, vrtwrk,
171  & vrtwrk(i1), vrtwrk(i2), vrtwrk(i3), vrtwrk(i4),
172  & vrtwrk(i5), vrtwrk(i6), VALUE, error, funcls, inform )
173  ELSE
174  funcls = 0
175  ENDIF
176 *
177 ****END SMPINT
178 *
179  END
180  SUBROUTINE smpchc( ND, NF, MINCLS, MAXCLS, EPSABS, EPSREL, SBRGNS,
181  & KEY, WRKLEN, RESTAR, RULCLS, MAXSUB, INFORM )
182 *
183 ****BEGIN PROLOGUE SMPCHC
184 ****AUTHOR
185 *
186 * Alan Genz
187 * Department of Mathematics
188 * Washington State University
189 * Pullman, WA 99164-3113, USA
190 *
191 ****LAST MODIFICATION 2001-07
192 ****PURPOSE SMPCHC checks validity of input parameters for SMPINT.
193 ****DESCRIPTION
194 * SMPCHC computes MAXSUB, RULCLS and INFORM as functions of
195 * input parameters for SMPINT, and checks the validity of
196 * input parameters for SMPINT.
197 *
198 * ON ENTRY
199 *
200 * ND Integer, number of variables, ND > 1.
201 * NF Integer, number of components of the integral.
202 * MINCLS Integer, minimum number of new FUNSUB calls.
203 * MAXCLS Integer, maximum number of new FUNSUB calls.
204 * The number of function values for each subregion is RULCLS.
205 * If
206 * KEY = 0, RULCLS = (ND+4)*(ND+3)*(ND+2)/6 + (ND+2)*(ND+1);
207 * KEY = 1, RULCLS = 2*ND+3;
208 * KEY = 2, RULCLS = (ND+3)*(ND+2)/2 + 2*(ND+1);
209 * KEY = 3, RULCLS = (ND+4)*(ND+3)*(ND+2)/6 + (ND+2)*(ND+1);
210 * KEY = 4, RULCLS = (ND+5)*(ND+4)*(ND+3)*(ND+2)/24
211 * + 5*(ND+2)*(ND+1)/2 .
212 * DIFCLS = 1 + 2*ND*( ND + 1 ).
213 * If RESTAR = 0, MAXCLS must be >= MAX(SBRGNS*RULCLS,MINCLS).
214 * If RESTAR = 1, MAXCLS must be >= MAX(4*RULCLS+DIFCLS,MINCLS).
215 * EPSABS Real, requested absolute accuracy.
216 * EPSREL Real, requested relative accuracy.
217 * SBRGNS Integer, initial number of simplices.
218 * KEY Integer, key to selected local integration rule.
219 * KEY = 0 gives the user a (default)degree 7 integration rule.
220 * KEY = 1 gives the user a degree 3 integration rule.
221 * KEY = 2 gives the user a degree 5 integration rule.
222 * KEY = 3 gives the user a degree 7 integration rule.
223 * KEY = 4 gives the user a degree 9 integration rule.
224 * WRKLEN Integer, length of the working array WORK.
225 * WRKLEN should be >= WRKSBS*( ND*(ND+1) + 2*NF + 3 )
226 * + (ND+1)*(ND+2) + 7*NF, where
227 * WRKSBS = SBRGNS + 3*( MAXCLS/RULCLS - SBRGNS*(1-RESTAR) )/4.
228 * RESTAR Integer.
229 * If RESTAR = 0, this is the first attempt to compute
230 * the integral over the SBRGNS input simplices.
231 * If RESTAR = 1, then we restart a previous attempt.
232 *
233 * ON RETURN
234 *
235 * RULCLS Integer, number of function values for each subregion.
236 * If
237 * KEY = 0, RULCLS = (ND+4)*(ND+3)*(ND+2)/6 + (ND+2)*(ND+1);
238 * KEY = 1, RULCLS = 2*ND+3;
239 * KEY = 2, RULCLS = (ND+3)*(ND+2)/2 + 2*(ND+1);
240 * KEY = 3, RULCLS = (ND+4)*(ND+3)*(ND+2)/6 + (ND+2)*(ND+1);
241 * KEY = 4, RULCLS = (ND+5)*(ND+4)*(ND+3)*(ND+2)/24
242 * + 5*(ND+2)*(ND+1)/2 .
243 * MAXSUB Integer, maximum allowed number of subregions for the
244 * given values of MAXCLS, WRKLEN, KEY and ND.
245 * INFORM Integer.
246 * INFORM = 0 for normal exit,
247 * INFORM = 2 if KEY < 0 or KEY > 4,
248 * INFORM = 3 if ND < 2,
249 * INFORM = 4 if NF < 1,
250 * INFORM = 5 if EPSABS < 0 and EPSREL < 0.,
251 * INFORM = 6 if WRKLEN is too small,
252 * INFORM = 7 if RESTAR < 0 or RESTAR > 1,
253 * INFORM = 8 if SBRGNS <= 0.
254 *
255 ****END PROLOGUE SMPCHC
256 *
257 * Global variables.
258 *
259  INTEGER ND, NF, MINCLS, MAXCLS, KEY, MAXSUB
260  INTEGER WRKLEN, INFORM, RESTAR, RULCLS, SBRGNS
261  DOUBLE PRECISION EPSABS, EPSREL
262 *
263 * Local variables.
264 *
265  INTEGER WRKDIF, DIFCLS
266 *
267 ****FIRST PROCESSING STATEMENT SMPCHC
268 *
269  inform = 0
270 *
271 * Check valid KEY.
272 *
273  IF ( key .LT. 0 .OR. key .GT. 4 ) inform = 2
274 *
275 * Check valid ND.
276 *
277  IF ( nd .LT. 2 ) inform = 3
278 *
279 * Check positive NF.
280 *
281  IF ( nf .LT. 1 ) inform = 4
282 *
283 * Check valid accuracy requests.
284 *
285  IF ( epsabs .LT. 0 .AND. epsrel .LT. 0 ) inform = 5
286 *
287 * Check workspace.
288 *
289  wrkdif = (nd+1)*(nd+2) + 7*nf
290  maxsub = ( wrklen - wrkdif )/( (nd+1)*nd + 2*nf + 3 )
291  IF ( maxsub .LE. sbrgns ) inform = 6
292 *
293 * Check valid RESTAR.
294 *
295  IF ( restar .NE. 0 .AND. restar .NE. 1 ) inform = 7
296 *
297 * Check valid SBRGNS.
298 *
299  IF ( sbrgns .LE. 0 ) inform = 8
300 *
301 * Compute RULCLS as a function of KEY and ND and check MAXCLS.
302 *
303  IF ( inform .EQ. 0 ) THEN
304  difcls = 1 + 2*nd*( nd + 1 )
305  IF (key .EQ. 0) rulcls = (nd+4)*(nd+3)*(nd+2)/6 + (nd+2)*(nd+1)
306  IF (key .EQ. 1) rulcls = 2*nd + 3
307  IF (key .EQ. 2) rulcls = (nd+3)*(nd+2)/2 + 2*(nd+1)
308  IF (key .EQ. 3) rulcls = (nd+4)*(nd+3)*(nd+2)/6 + (nd+2)*(nd+1)
309  IF (key .EQ. 4) rulcls = (nd+5)*(nd+4)*(nd+3)*(nd+2)/24
310  & + 5*(nd+2)*(nd+1)/2
311  IF ( restar.EQ.0 .AND. maxcls.LT.max(sbrgns*rulcls,mincls) .OR.
312  & restar.EQ.1 .AND. maxcls.LT.max(4*rulcls+difcls,mincls) )
313  & inform = 1
314  ENDIF
315 *
316 ****END SMPCHC
317 *
318  END
319 *
320  SUBROUTINE smpsad( ND, NF, FUNSUB, MINCLS, MAXCLS, EPSABS, EPSREL,
321  & RESTAR, KEY, RULCLS, MAXSUB, SBRGNS, VERTCS, VALUES, ERRORS,
322  & VOLUMS, GREATS, PONTRS, WORK, VALUE, ERROR, FUNCLS, INFORM )
323 *
324 ****BEGIN PROLOGUE SMPSAD
325 ****KEYWORDS automatic multidimensional integrator,
326 * n-dimensional simplex,
327 * general purpose, global adaptive
328 ****AUTHOR
329 *
330 * Alan Genz
331 * Department of Mathematics
332 * Washington State University
333 * Pullman, WA 99164-3113, USA
334 *
335 ****LAST MODIFICATION 2001-07
336 ****PURPOSE The routine calculates an approximation to a given
337 * vector of definite integrals, I, over a hyper-rectangular
338 * region hopefully satisfying for each component of I the
339 * following claim for accuracy:
340 * ABS( I(K) - VALUE(K) ) .LE. MAX( EPSABS, EPSREL*ABS(I(K) ) )
341 ****DESCRIPTION Computation of integrals over hyper-rectangular regions.
342 * SMPSAD repeatedly subdivides the regions of integration
343 * and estimates the integrals and the errors over the
344 * subregions with greatest estimated errors until the error
345 * request is met or MAXSUB subregions are stored. The regions
346 * are divided into three or four equally sized parts along
347 * the direction(s) with greatest absolute fourth difference.
348 *
349 * ON ENTRY
350 *
351 * ND Integer, number of variables, ND > 1.
352 * NF Integer, number of components of the integral.
353 * FUNSUB Externally declared subroutine for computing components of
354 * the integrand at the given evaluation point.
355 * It must have parameters (ND,X,NF,FUNVLS)
356 * Input parameters:
357 * ND Integer that gives the dimension of I
358 * X Real array of dimension ND that contains the
359 * evaluation point.
360 * NF Integer that gives the number of components of I.
361 * Output parameter:
362 * FUNVLS Real array of dimension NF that contains the
363 * components of the integrand.
364 * MINCLS Integer.
365 * The computations proceed until there are at least
366 * MINCLS FUNSUB calls.
367 * MAXCLS Integer.
368 * The computations proceed until further subdivision would
369 * require more than MAXCLS FUNSUB calls. When RESTAR = 1,
370 * this is the number of new FUNSUB calls.
371 * EPSABS Real, requested absolute accuracy.
372 * EPSREL Real, requested relative accuracy.
373 * RESTAR Integer.
374 * If RESTAR = 0, this is the first attempt to compute
375 * the integral.
376 * If RESTAR = 1, then we restart a previous attempt.
377 * (In this case the output parameters from SMPSAD
378 * must not be changed since the last exit from SMPSAD.)
379 * KEY Integer, key to selected local integration rule.
380 * KEY = 1 gives the user a degree 3 integration rule.
381 * KEY = 2 gives the user a degree 5 integration rule.
382 * KEY = 3 gives the user a degree 7 integration rule.
383 * KEY = 4 gives the user a degree 9 integration rule.
384 * RULCLS Integer, number of FUNSUB calls needed for each subregion.
385 * MAXSUB Integer; computations proceed until there are at most
386 * MAXSUB subregions in the data structure.
387 * SBRGNS Integer.
388 * If RESTAR = 0, then SBRGNS must specify the number
389 * of subregions stored in a previous call to SMPSAD.
390 * VERTCS Real array of dimension (ND,0:ND,*).
391 * Simplex vertices for each subregion; for subregion K vertex
392 * J must have components VERTEX(I,J,K), I = 1, 2, ..., ND.
393 * VALUES Real array of dimension (NF,*), for estimated values of the
394 * integrals over the subregions.
395 * ERRORS Real array of dimension (NF,*).
396 * Used to store the corresponding estimated errors.
397 * Used to store the half widths of the stored subregions.
398 * GREATS Real array of dimension (*).
399 * Used to store the greatest estimated errors in subregions.
400 * PONTRS Real array of dimension (*), for the pointers from the
401 * subregion heap to the actual subregions.
402 * WORK Real array, used in SMPVOL, SMPDFS, and SMPRUL.
403 *
404 * ON RETURN
405 *
406 * SBRGNS Integer, number of stored subregions.
407 * VALUE Real array of dimension NF.
408 * Approximations to all components of the integral.
409 * ERROR Real array of dimension NF.
410 * Estimates of absolute accuracies.
411 * FUNCLS Integer, number of new FUNSUB calls used by SMPSAD.
412 * INFORM Integer.
413 * INFORM = 0 for normal exit, when ERROR(K) <= EPSABS or
414 * ERROR(K) <= ABS(VALUE(K))*EPSREL, 1 <= K <= NF,
415 * with MAXSUB or fewer subregions processed.
416 * INFORM = 1 if MAXSUB was too small for SMPSAD
417 * to obtain the required accuracy. In this case SMPSAD
418 * returns values of VALUE with estimated absolute
419 * accuracies ERROR.
420 *
421 ****REFERENCES
422 ****ROUTINES CALLED SMPSTR, SMPVOL, SMPRUL
423 ****END PROLOGUE SMPSAD
424 *
425 * Global variables.
426 *
427  EXTERNAL funsub
428  INTEGER ND, NF, RULCLS, MINCLS, MAXCLS, MAXSUB, KEY, RESTAR
429  INTEGER FUNCLS, SBRGNS, INFORM
430  DOUBLE PRECISION EPSABS, EPSREL, VALUE(NF), ERROR(NF)
431  DOUBLE PRECISION VALUES(NF,*), ERRORS(NF,*), VERTCS(ND,0:ND,*)
432  DOUBLE PRECISION VOLUMS(*), GREATS(*), PONTRS(*), WORK(*)
433 *
434 * Local variables.
435 *
436 *
437 * MXNWSB is the maxiumum number of new subregions per subdivision.
438 *
439  INTEGER I, INDEX, J, TOP, MXNWSB, NEWSBS, DFCOST, RGNCLS
440  PARAMETER ( MXNWSB = 4 )
441  double precision smpvol, tune
442  parameter( tune = 1 )
443 *
444 ****FIRST PROCESSING STATEMENT SMPSAD
445 *
446 *
447 * Initialize for rule parameters.
448 *
449  funcls = 0
450  dfcost = 1 + 2*nd*( nd + 1 )
451 *
452 * If RESTAR = 0, initialize for first call.
453 *
454  IF ( restar .EQ. 0 ) THEN
455 *
456 * Initialize FUNCLS, and VALUE and ERROR arrays.
457 *
458  DO j = 1, nf
459  value(j) = 0
460  error(j) = 0
461  END DO
462  DO index = 1, sbrgns
463 *
464 * Call SMPVOL to compute the simplex volume(s).
465 *
466  volums(index) = smpvol( nd, vertcs(1,0,index), work )
467 *
468 * Apply basic rule over each simplex.
469 *
470  CALL smprul( tune, nd, vertcs(1,0,index), volums(index),
471  & nf, funsub, key, values(1,index), errors(1,index),
472  & greats(index), work, work(2*nd+2) )
473 *
474 * Add new contributions to VALUE and ERROR.
475 * Store results in heap.
476 *
477  DO j = 1, nf
478  value(j) = value(j) + values(j,index)
479  error(j) = error(j) + errors(j,index)
480  END DO
481  CALL smpstr( index, index, pontrs, greats )
482  funcls = funcls + rulcls
483  END DO
484  END IF
485  inform = max( 0, min( mincls - funcls, 1 ) )
486  DO j = 1, nf
487  IF( error(j) .GT. max(epsabs,epsrel*abs(value(j))) ) inform = 1
488  END DO
489 *
490 * End initialisation.
491 *
492  DO WHILE ( inform .GT. 0 .AND. sbrgns + mxnwsb - 1 .LE. maxsub
493  & .AND. funcls + dfcost + mxnwsb*rulcls .LE. maxcls )
494 *
495 * Begin loop while error is too large, and FUNCLS and SBRGNS
496 * are not too large.
497 *
498 * Adjust VALUE and ERROR.
499 *
500  top = pontrs(1)
501  DO j = 1, nf
502  value(j) = value(j) - values(j,top)
503  error(j) = error(j) - errors(j,top)
504  END DO
505 *
506 * Determine NEWSBS new subregions.
507 *
508  CALL smpdfs( nd, nf, funsub, top, sbrgns, vertcs,
509  & volums, work, work(nd+1), work(2*nd+1),
510  & work(3*nd+1), work(3*nd+1+5*nf), newsbs )
511 *
512 * Apply basic rule, store results in heap and
513 * add new contributions to VALUE and ERROR.
514 *
515  index = top
516  DO i = 1, newsbs
517  CALL smprul( tune, nd, vertcs(1,0,index), volums(index), nf,
518  & funsub, key, values(1,index), errors(1,index),
519  & greats(index), work, work(2*nd+2) )
520  CALL smpstr( index, sbrgns+i-1, pontrs, greats )
521  DO j = 1, nf
522  value(j) = value(j) + values(j,index)
523  error(j) = error(j) + errors(j,index)
524  END DO
525  index = sbrgns + i
526  END DO
527  funcls = funcls + dfcost + newsbs*rulcls
528  sbrgns = sbrgns + newsbs - 1
529 *
530 * Check for error termination.
531 *
532  inform = max( 0, min( mincls - funcls, 1 ) )
533  DO j = 1, nf
534  IF( error(j) .GT. max(epsabs,epsrel*abs(value(j))) )
535  & inform = 1
536  END DO
537  END DO
538 *
539 * Compute more accurate values of VALUE and ERROR.
540 *
541  DO i = 1, nf
542  value(i) = 0
543  error(i) = 0
544  DO j = 1, sbrgns
545  value(i) = value(i) + values(i,j)
546  error(i) = error(i) + errors(i,j)
547  END DO
548  END DO
549 *
550 ****END SMPSAD
551 *
552  END
553 *
554  DOUBLE PRECISION FUNCTION smpvol( ND, VERTEX, WORK )
555 *
556 ****BEGIN PROLOGUE SMPVOL
557 ****KEYWORDS simplex volume
558 ****PURPOSE Function to compute the scaled volume of a simplex.
559 ****AUTHOR
560 *
561 * Alan Genz
562 * Department of Mathematics
563 * Washington State University
564 * Pullman, WA 99164-3113, USA
565 *
566 ****LAST MODIFICATION 96-12
567 ****DESCRIPTION SMPVOL computes the volume of an ND-simplex scaled by
568 * using Gauss elimination to compute a determinant.
569 *
570 * ON ENTRY
571 *
572 * ND Integer, number of variables.
573 * VERTEX Real array of dimension (ND,0:ND), of simplex vertices;
574 * vertex J must have components VERTEX(I,J), I = 1, 2, ..., ND.
575 * WORK Real work array of dimension at least ND*ND.
576 *
577 * ON RETURN
578 *
579 * SMPVOL Real, value for the volume.
580 *
581 ****ROUTINES CALLED: None.
582 *
583 ****END PROLOGUE SMPVOL
584 *
585 * Global variables.
586 *
587  INTEGER ND
588  DOUBLE PRECISION VERTEX(ND,0:ND), WORK(ND,*)
589 *
590 * Local variables.
591 *
592  INTEGER I, J, K, PIVPOS
593  DOUBLE PRECISION MULT, VOL, WTEMP
594 *
595 ****FIRST PROCESSING STATEMENT SMPVOL
596 *
597 *
598 * Copy vertex differences to WORK array.
599 *
600  DO j = 1, nd
601  DO i = 1, nd
602  work(i,j) = vertex(i,j) - vertex(i,0)
603  END DO
604  END DO
605 *
606 * Use Gauss elimination with partial pivoting.
607 *
608  vol = 1
609  DO k = 1, nd
610  pivpos = k
611  DO j = k+1, nd
612  IF ( abs(work(k,j)) .GT. abs(work(k,pivpos)) ) pivpos = j
613  END DO
614  DO i = k, nd
615  wtemp = work(i,k)
616  work(i,k) = work(i,pivpos)
617  work(i,pivpos) = wtemp
618  END DO
619  vol = vol*work(k,k)/k
620  DO j = k+1, nd
621  mult = work(k,j)/work(k,k)
622  DO i = k+1, nd
623  work(i,j) = work(i,j) - mult*work(i,k)
624  END DO
625  END DO
626  END DO
627  smpvol = abs(vol)
628 *
629 ****END SMPVOL
630 *
631  END
632 *
633  SUBROUTINE smprul( TUNE, ND, VERTEX, VOLUME, NF, INTGND,
634  & INKEY, BASVAL, RGNERR, GREAT, GT, RULE )
635 *
636 ****BEGIN PROLOGUE SMPRUL
637 ****KEYWORDS basic numerical integration rule
638 ****PURPOSE To compute basic integration rule values.
639 ****AUTHOR
640 *
641 * Alan Genz
642 * Department of Mathematics
643 * Washington State University
644 * Pullman, WA 99164-3113, USA
645 * AlanGenz@wsu.edu
646 *
647 ****LAST MODIFICATION 2001-07
648 ****DESCRIPTION SMPRUL computes basic integration rule values for a
649 * vector of integrands over a hyper-rectangular region.
650 * These are estimates for the integrals. SMPRUL also computes
651 * estimates for the errors.
652 *
653 * ON ENTRY
654 *
655 * TUNE Real, tuning parameter, with 0 <= TUNE <= 1, with
656 * TUNE = 1 for the most conservative error estimates.
657 * If TUNE < 0, only the rule parameters are computed.
658 * ND Integer, number of variables.
659 * VERTEX Real array of dimension (ND,0:ND).
660 * The simplex vertices; vertex J must have components
661 * VERTEX(I,J), I = 1, 2, ..., ND.
662 * NF Integer, number of components for the vector integrand.
663 * INTGND Subroutine for computing components of the integrand at Z.
664 * It must have parameters (ND,X,NF,FUNVLS)
665 * Input parameters:
666 * ND Integer that gives the dimension.
667 * X Real array of dimension ND that contains the
668 * evaluation point.
669 * NF Integer that gives the number of components of I.
670 * Output parameter:
671 * FUNVLS Real array of dimension NF that contains the
672 * components of the integrand.
673 * INKEY Integer rule parameter.
674 * If INKEY .GT. 0 and INKEY .LT. 5 then a rule of degree
675 * 2*INKEY + 1; otherwise default degree 7 rule is used.
676 * GT Real work array of length 2*ND+1.
677 * RULE Real work array of dimension (NF,7).
678 *
679 * ON RETURN
680 *
681 * BASVAL Real array of length NF, values for the basic rule for
682 * each component of the integrand.
683 * RGNERR Real array of length NF, error estimates for BASVAL.
684 * GREAT Real, maximum component of RGNERR.
685 *
686 *
687 ****ROUTINES CALLED: SMPRMS, SYMRUL
688 *
689 ****END PROLOGUE SMPRUL
690 *
691 * Global variables.
692 *
693  EXTERNAL intgnd
694  INTEGER NF, ND, INKEY
695  DOUBLE PRECISION VERTEX(ND,0:ND), BASVAL(NF), RGNERR(NF)
696  DOUBLE PRECISION VOLUME, TUNE, GREAT
697 *
698 * Local variables.
699 *
700 * WTS Integer number of weights in the integration rules.
701 * W Real array of dimension (WTS,RLS).
702 * The weights for the basic and null rules.
703 * W(1,1),...,W(WTS,1) are weights for the basic rule.
704 * W(1,I),...,W(WTS,I), for I > 1 are null rule weights.
705 * G Real array of dimension (0:4, WTS).
706 * The fully symmetric sum generators for the rules.
707 *
708  INTEGER KEY, NUMNUL, RLS, WTS, MXW, MXRLS, MXG
709  PARAMETER( MXW = 21, mxrls = 7, mxg = 4 )
710  DOUBLE PRECISION W( MXW, MXRLS ), G( 0:MXG, MXW ), WTSUM
711  DOUBLE PRECISION GT( 0:2*ND ), RULE( NF, MXRLS )
712  DOUBLE PRECISION NORMCF, NORMNL, NORMCP, ALPHA(MXRLS)
713  DOUBLE PRECISION RATIO, ERRCOF, RATMIN, SMALL, SMPROD
714  parameter( ratmin = 1d-1, small = 1d-12 )
715  INTEGER I, J, K, OLDKEY, OLDN, PTS(MXW)
716  SAVE oldkey, oldn, key, pts, w, g, rls, wts
717  DATA oldkey, oldn/ -1, 0 /
718 *
719 ****FIRST PROCESSING STATEMENT SMPRUL
720 *
721  IF ( oldkey .NE. inkey .OR. oldn .NE. nd ) THEN
722  oldn = nd
723  oldkey = inkey
724  IF ( inkey .GT. 0 .AND. inkey .LT. 5 ) THEN
725  key = inkey
726  ELSE
727  key = 3
728  END IF
729 *
730 * Compute WTS, RLS, weights, generators, ERRCOF and PTS.
731 *
732  CALL smprms( nd, key, mxw, w, mxg, g, wts, rls, pts )
733 *
734 * Orthogonalize and normalize null rules.
735 *
736  normcf = smprod( wts, pts, w(1,1), w(1,1) )
737  DO k = 2, rls
738  DO j = 2, k-1
739  alpha(j) = -smprod( wts, pts, w(1,j), w(1,k) )
740  END DO
741  DO i = 1, wts
742  wtsum = 0
743  DO j = 2, k-1
744  wtsum = wtsum + w(i,j)*alpha(j)
745  END DO
746  w(i,k) = w(i,k) + wtsum/normcf
747  END DO
748  normnl = smprod( wts, pts, w(1,k), w(1,k) )
749  DO i = 1, wts
750  w(i,k) = w(i,k)*sqrt( normcf/normnl )
751  END DO
752  END DO
753  ENDIF
754  IF ( tune .GE. 0 ) THEN
755 *
756 * Compute the rule values.
757 *
758  DO i = 1, nf
759  DO j = 1, rls
760  rule(i,j) = 0
761  END DO
762  END DO
763  DO k = 1, wts
764  IF ( pts(k) .GT. 0 ) THEN
765  DO i = 0, min(nd,mxg-1)
766  gt(i) = g(i,k)
767  END DO
768  IF ( nd .GE. mxg ) CALL smpcpy( mxg, nd, gt, g(mxg,k) )
769  CALL smpsms( nd, vertex, nf, intgnd, gt, basval,
770  & gt(nd+1), rgnerr )
771  DO j = 1, rls
772  DO i = 1, nf
773  rule(i,j) = rule(i,j) + w(k,j)*basval(i)
774  END DO
775  END DO
776  END IF
777  END DO
778 *
779 * Scale integral values and compute the error estimates.
780 *
781  errcof = ( 8*tune + ( 1 - tune ) )
782  great = 0
783  DO i = 1, nf
784  basval(i) = rule(i,1)
785  normcf = abs( basval(i) )
786  rgnerr(i) = 0
787  ratio = ratmin
788  DO k = rls, 3, -2
789  normnl = max( abs( rule(i,k) ), abs( rule(i,k-1) ) )
790  IF ( normnl .GT. small*normcf .AND. k .LT. rls )
791  & ratio = max( normnl/normcp, ratio )
792  rgnerr(i) = max( normnl, rgnerr(i) )
793  normcp = normnl
794  END DO
795  IF( ratio .GE. 1 ) THEN
796  rgnerr(i) = tune*rgnerr(i) + ( 1 - tune )*normcp
797  ELSE IF ( key .GT. 1 ) THEN
798  rgnerr(i) = ratio*normcp
799  END IF
800  rgnerr(i) = volume*max( errcof*rgnerr(i), small*normcf )
801  basval(i) = volume*basval(i)
802  great = max( great, rgnerr(i) )
803  END DO
804  END IF
805 *
806 ****END SMPRUL
807 *
808  END
809 *
810  DOUBLE PRECISION FUNCTION smprod( N, W, X, Y )
811  INTEGER N, I, W(*)
812  DOUBLE PRECISION X(*), Y(*), SUM
813  SUM = 0
814  DO i = 1, n
815  sum = sum + w(i)*x(i)*y(i)
816  END DO
817  smprod = sum
818  END
819 *
820  SUBROUTINE smprms( ND, KEY, MXW, W, MXG, G, WTS, RLS, PTS )
821 *
822 ****BEGIN PROLOGUE SMPRMS
823 ****KEYWORDS basic integration rule, degree 2*KEY+1
824 ****PURPOSE To initialize a degree 2*KEY+1 basic rule and null rules.
825 ****AUTHOR
826 *
827 * Alan Genz
828 * Department of Mathematics
829 * Washington State University
830 * Pullman, WA 99164-3113, USA
831 * AlanGenz@wsu.edu
832 *
833 ****LAST MODIFICATION 2001-07
834 ****DESCRIPTION SMPRMS initializes a degree 2*KEY+1 rule, and
835 * and max(2*KEY,2) lower degree null rules.
836 *
837 * ON ENTRY
838 *
839 * ND Integer, number of variables.
840 * KEY Integer, < 5 and >= 0, rule parameter.
841 * If KEY > 0 a degree 2*KEY+1 rule is initialized.
842 * If KEY = 0 a degree 7 rule is initialized.
843 *
844 * ON RETURN
845 * RLS Integer, total number of rules.
846 * WTS Integer, total number of weights in each of the rules.
847 * W Real array of dimension (MXW,*).
848 * The weights for the basic and null rules.
849 * W(1,1),...,W(WTS,1) are weights for the basic rule.
850 * W(I,1),...,W(WTS,I) for I .GT. 1 are null rule weights.
851 * G Real array of dimension (0:MXG,MXW).
852 * The fully symmetric sum generators for the rules.
853 * G(0,J), ..., G(MXG,J) are the generators for the
854 * points associated with the Jth weights.
855 * PTS Integer array of length (MXW). PTS(J) is number of integrand
856 * values needed for generator J.
857 *
858 ****REFERENCES
859 *
860 * Axel Grundmann and H. M. Moller
861 * "Invariant Integration Formulas for the n-Simplex by Combinatorial
862 * Methods", SIAM J Numer. Anal. 15(1978), 282--290,
863 * and
864 * A. H. Stroud
865 * "A Fifth Degree Integration Formula for the n-Simplex
866 * SIAM J Numer. Anal. 6(1969), 90--98,
867 * and
868 * I. P. Mysovskikh
869 * "On a cubature formula for the simplex"
870 * Vopros. Vycisl. i Prikl. Mat., Tashkent 51(1978), 74--90.
871 *
872 *
873 ****ROUTINES CALLED NONE
874 ****END PROLOGUE SMPRMS
875 *
876 * Global variables
877 *
878  INTEGER ND, KEY, WTS, MXW, RLS, MXG
879  INTEGER PTS(MXW)
880  DOUBLE PRECISION W(MXW,*), G(0:MXG,*)
881 *
882 * Local Variables
883 *
884  DOUBLE PRECISION ONE, FFTEEN
885  parameter( one = 1, ffteen = 15 )
886  DOUBLE PRECISION DR, DR2, DR4, DR6, DR8
887  DOUBLE PRECISION R1, S1, R2, S2, U1, V1, U2, V2, L1, L2, D1, D2
888  DOUBLE PRECISION A1, A2, A3, P0, P1, P2, P3, U5, U6, U7, SG
889  DOUBLE PRECISION R, A, P, Q, TH, TP
890  INTEGER IW, GMS, I, J
891 *
892 ****FIRST PROCESSING STATEMENT SMPRMS
893 *
894 *
895 * Initialize RLS and GMS.
896 *
897  IF ( key .EQ. 1 ) THEN
898  rls = 3
899  gms = 2
900  wts = 3
901  ELSE IF ( key .EQ. 2 ) THEN
902  rls = 5
903  gms = 4
904  wts = 6
905  ELSE IF ( key .EQ. 3 .OR. key .EQ. 0 ) THEN
906  rls = 7
907  gms = 7
908  wts = 11
909  ELSE IF ( key .EQ. 4 ) THEN
910  rls = 7
911  IF ( nd .EQ. 2 ) THEN
912  gms = 11
913  wts = 20
914  ELSE
915  gms = 12
916  wts = 21
917  END IF
918  END IF
919 *
920 * Initialize generators, weights and PTS.
921 *
922  DO i = 1, wts
923  DO j = 1, rls
924  w(i,j) = 0
925  END DO
926  pts(i) = 0
927  END DO
928 *
929 * Compute generator, PTS and weight values for all rules.
930 *
931  dr = nd
932  dr2 = ( dr + 1 )*( dr + 2 )
933  dr4 = dr2*( dr + 3 )*( dr + 4 )
934  dr6 = dr4*( dr + 5 )*( dr + 6 )
935  dr8 = dr6*( dr + 7 )*( dr + 8 )
936  CALL smpcpy( 0, mxg, g(0,1), 1/( dr + 1 ) )
937  pts(1) = 1
938  r1 = ( dr + 4 - sqrt(ffteen) )/( dr*dr + 8*dr + 1 )
939  s1 = 1 - dr*r1
940  l1 = s1 - r1
941  g(0 ,gms+1) = s1
942  CALL smpcpy( 1, mxg, g(0,gms+1), r1 )
943  DO i = 1, mxg
944  g(i,gms+1) = r1
945  END DO
946  pts(gms+1) = dr + 1
947  iw = rls
948  IF ( key .LT. 4 ) THEN
949 *
950 * Compute weights for special degree 1 rule.
951 *
952  w(1,iw) = 1
953  iw = iw - 1
954  w(gms+1,iw) = 1/( dr + 1 )
955  iw = iw - 1
956  END IF
957 *
958 * Compute weights, generators and PTS for degree 3 rule.
959 *
960  g(0,2) = 3/( dr + 3 )
961  CALL smpcpy( 1, mxg, g(0,2), 1/( dr + 3 ) )
962  pts(2) = dr + 1
963  w(2,iw) = ( dr + 3 )**3/( 4*dr2*( dr + 3 ) )
964  IF ( key .GT. 1 ) THEN
965  iw = iw - 1
966 *
967 * Compute weights, generators and PTS for degree 3 and 5 rules.
968 *
969  IF ( nd .EQ. 2 ) THEN
970 *
971 * Special degree 3 rule.
972 *
973  l2 = .62054648267200632589046034361711d0
974  l1 = -sqrt( one/2 - l2**2 )
975  r1 = ( 1 - l1 )/3
976  s1 = 1 - 2*r1
977  g(0,gms+1) = s1
978  CALL smpcpy( 1, mxg, g(0,gms+1), r1 )
979  pts(gms+1) = 3
980  w(gms+1,iw) = one/6
981  r2 = ( 1 - l2 )/3
982  s2 = 1 - 2*r2
983  g(0,gms+2) = s2
984  CALL smpcpy( 1, mxg, g(0,gms+2), r2 )
985  pts(gms+2) = 3
986  w(gms+2,iw) = one/6
987  ELSE
988 *
989 * Degree 3 rule using Stroud points.
990 *
991  r2 = ( dr + 4 + sqrt(ffteen) )/( dr*dr + 8*dr + 1 )
992  s2 = 1 - dr*r2
993  l2 = s2 - r2
994  g(0,gms+2) = s2
995  CALL smpcpy( 1, mxg, g(0,gms+2), r2 )
996  pts(gms+2) = dr + 1
997  w(gms+2,iw) = ( 2/(dr+3) - l1 )/( dr2*(l2-l1)*l2**2 )
998  w(gms+1,iw) = ( 2/(dr+3) - l2 )/( dr2*(l1-l2)*l1**2 )
999  END IF
1000  iw = iw - 1
1001 *
1002 * Grundmann-Moller degree 5 rule.
1003 *
1004  g(0,3) = 5/( dr + 5 )
1005  CALL smpcpy( 1, mxg, g(0,3), 1/( dr + 5 ) )
1006  pts(3) = dr + 1
1007  g(0,4) = 3/( dr + 5 )
1008  g(1,4) = 3/( dr + 5 )
1009  CALL smpcpy( 2, mxg, g(0,4), 1/( dr + 5 ) )
1010  pts(4) = ( ( dr + 1 )*dr )/2
1011  w(2,iw) = -( dr + 3 )**5/( 16*dr4 )
1012  w(3,iw) = ( dr + 5 )**5/( 16*dr4*( dr + 5 ) )
1013  w(4,iw) = ( dr + 5 )**5/( 16*dr4*( dr + 5 ) )
1014  END IF
1015  IF ( key .GT. 2 ) THEN
1016  iw = iw - 1
1017 *
1018 * Compute weights, generators and PTS for degree 5 and 7 rules.
1019 *
1020 *
1021 * Stroud degree 5 rule.
1022 *
1023  u1 = ( dr + 7 + 2*sqrt(ffteen) )/( dr*dr + 14*dr - 11 )
1024  v1 = ( 1 - ( dr - 1 )*u1 )/2
1025  d1 = v1 - u1
1026  g(0,gms+3) = v1
1027  g(1,gms+3) = v1
1028  CALL smpcpy( 2, mxg, g(0,gms+3), u1 )
1029  pts(gms+3) = ( ( dr + 1 )*dr )/2
1030  u2 = ( dr + 7 - 2*sqrt(ffteen) )/( dr*dr + 14*dr - 11 )
1031  v2 = ( 1 - ( dr - 1 )*u2 )/2
1032  d2 = v2 - u2
1033  g(0,gms+4) = v2
1034  g(1,gms+4) = v2
1035  CALL smpcpy( 2, mxg, g(0,gms+4), u2 )
1036  pts(gms+4) = ( ( dr + 1 )*dr )/2
1037  IF ( nd .EQ. 2 ) THEN
1038  w(gms+3,iw) = ( 155 - sqrt(ffteen) )/1200
1039  w(gms+4,iw) = ( 155 + sqrt(ffteen) )/1200
1040  w(1, iw) = 1 - 3*( w(gms+3,iw) + w(gms+4,iw) )
1041  ELSE IF ( nd .EQ. 3 ) THEN
1042  w(gms+1,iw) = ( 2665 + 14*sqrt(ffteen) )/37800
1043  w(gms+2,iw) = ( 2665 - 14*sqrt(ffteen) )/37800
1044  w(gms+3,iw) = 2*ffteen/567
1045  pts(gms+4) = 0
1046  ELSE
1047  w(gms+1,iw) = ( 2*(27-dr)/(dr+5)-l2*(13-dr) )
1048  & /( l1**4*(l1-l2)*dr4 )
1049  w(gms+2,iw) = ( 2*(27-dr)/(dr+5)-l1*(13-dr) )
1050  & /( l2**4*(l2-l1)*dr4 )
1051  w(gms+3,iw)=( 2/( dr + 5 ) - d2 )/( dr4*( d1 - d2 )*d1**4 )
1052  w(gms+4,iw)=( 2/( dr + 5 ) - d1 )/( dr4*( d2 - d1 )*d2**4 )
1053  END IF
1054  iw = iw - 1
1055 *
1056 * Grundmann-Moller degree 7 rule.
1057 *
1058  g(0,5) = 7/( dr + 7 )
1059  CALL smpcpy( 1, mxg, g(0,5), 1/( dr + 7 ) )
1060  pts(5) = dr + 1
1061  g(0,6) = 5/( dr + 7 )
1062  g(1,6) = 3/( dr + 7 )
1063  CALL smpcpy( 2, mxg, g(0,6), 1/( dr + 7 ) )
1064  pts(6) = ( dr + 1 )*dr
1065  g(0,7) = 3/( dr + 7 )
1066  g(1,7) = 3/( dr + 7 )
1067  g(2,7) = 3/( dr + 7 )
1068  CALL smpcpy( 3, mxg, g(0,7), 1/( dr + 7 ) )
1069  pts(7) = ( ( dr + 1 )*dr*( dr - 1 ) )/6
1070  w(2,iw) = ( dr + 3 )**7/( 2*64*dr4*( dr + 5 ) )
1071  w(3,iw) = -( dr + 5 )**7/( 64*dr6 )
1072  w(4,iw) = -( dr + 5 )**7/( 64*dr6 )
1073  w(5,iw) = ( dr + 7 )**7/( 64*dr6*( dr + 7 ) )
1074  w(6,iw) = ( dr + 7 )**7/( 64*dr6*( dr + 7 ) )
1075  w(7,iw) = ( dr + 7 )**7/( 64*dr6*( dr + 7 ) )
1076  END IF
1077  IF ( key .EQ. 4 ) THEN
1078  iw = iw - 1
1079 *
1080 * Compute weights, generators and PTS for degree 7, 9 rules.
1081 *
1082 * Mysovskikh degree 7 rule.
1083 *
1084  sg = 1/( 23328*dr6 )
1085  u5 = -6**3*sg*( 52212 - dr*( 6353 + dr*( 1934 - dr*27 ) ) )
1086  u6 = 6**4*sg*( 7884 - dr*( 1541 - dr*9 ) )
1087  u7 = -6**5*sg*( 8292 - dr*( 1139 - dr*3 ) )/( dr + 7 )
1088  p0 = -144*( 142528 + dr*( 23073 - dr*115 ) )
1089  p1 = -12*( 6690556 + dr*( 2641189 + dr*( 245378 - dr*1495 ) ) )
1090  p2 = -16*(6503401 + dr*(4020794+dr*(787281+dr*(47323-dr*385))))
1091  p3 = -( 6386660 + dr*(4411997+dr*(951821+dr*(61659-dr*665))) )
1092  & *( dr + 7 )
1093  a = p2/( 3*p3 )
1094  p = a*( p1/p2 - a )
1095  q = a*( 2*a*a - p1/p3 ) + p0/p3
1096  r = sqrt( -p**3 )
1097  th = acos( -q/( 2*r ) )/3
1098  r = 2*r**( one/3 )
1099  tp = 2*acos(-one)/3
1100  a1 = -a + r*cos( th )
1101  a2 = -a + r*cos( th + tp + tp )
1102  a3 = -a + r*cos( th + tp )
1103  g(0,gms+5) = ( 1 - dr*a1 )/( dr + 1 )
1104  CALL smpcpy( 1, mxg, g(0,gms+5), ( 1 + a1 )/( dr + 1 ) )
1105  pts(gms+5) = dr + 1
1106  g(0,gms+6) = ( 1 - dr*a2 )/( dr + 1 )
1107  CALL smpcpy( 1, mxg, g(0,gms+6), ( 1 + a2 )/( dr + 1 ) )
1108  pts(gms+6) = dr + 1
1109  g(0,gms+7) = ( 1 - dr*a3 )/( dr + 1 )
1110  CALL smpcpy( 1, mxg, g(0,gms+7), ( 1 + a3 )/( dr + 1 ) )
1111  pts(gms+7) = dr + 1
1112  w(gms+5,iw) = ( u7-(a2+a3)*u6+a2*a3*u5 )
1113  & /( a1**2-(a2+a3)*a1+a2*a3 )/a1**5
1114  w(gms+6,iw) = ( u7-(a1+a3)*u6+a1*a3*u5 )
1115  & /( a2**2-(a1+a3)*a2+a1*a3 )/a2**5
1116  w(gms+7,iw) = ( u7-(a2+a1)*u6+a2*a1*u5 )
1117  & /( a3**2-(a2+a1)*a3+a2*a1 )/a3**5
1118  g(0,gms+8) = 4/( dr + 7 )
1119  g(1,gms+8) = 4/( dr + 7 )
1120  CALL smpcpy( 2, mxg, g(0,gms+8), 1/( dr + 7 ) )
1121  pts(gms+8) = ( ( dr + 1 )*dr )/2
1122  w(gms+8,iw) = 10*(dr+7)**6/( 729*dr6 )
1123  g(0,gms+9) = 11/( dr + 7 )/2
1124  g(1,gms+9) = 5/( dr + 7 )/2
1125  CALL smpcpy( 2, mxg, g(0,gms+9), 1/( dr + 7 ) )
1126  pts(gms+9) = ( ( dr + 1 )*dr )
1127  w(gms+9,iw) = 64*(dr+7)**6/( 6561*dr6 )
1128  w( 4,iw) = w(4,iw+1)
1129  w( 7,iw) = w(7,iw+1)
1130  iw = iw - 1
1131 *
1132 * Grundmann-Moller degree 9 rule.
1133 *
1134  g(0,8) = 9/( dr + 9 )
1135  CALL smpcpy( 1, mxg, g(0, 8), 1/( dr + 9 ) )
1136  pts(8) = dr + 1
1137  g(0,9) = 7/( dr + 9 )
1138  g(1,9) = 3/( dr + 9 )
1139  CALL smpcpy( 2, mxg, g(0, 9), 1/( dr + 9 ) )
1140  pts(9) = ( dr + 1 )*dr
1141  g(0,10) = 5/( dr + 9 )
1142  g(1,10) = 5/( dr + 9 )
1143  CALL smpcpy( 2, mxg, g(0,10), 1/( dr + 9 ) )
1144  pts(10) = ( ( dr + 1 )*dr )/2
1145  g(0,11) = 5/( dr + 9 )
1146  g(1,11) = 3/( dr + 9 )
1147  g(2,11) = 3/( dr + 9 )
1148  CALL smpcpy( 3, mxg, g(0,11), 1/( dr + 9 ) )
1149  pts(11) = ( ( dr + 1 )*dr*( dr - 1 ) )/2
1150  w(2 ,iw) = -( dr + 3 )**9/( 6*256*dr6 )
1151  w(3 ,iw) = ( dr + 5 )**9/( 2*256*dr6*( dr + 7 ) )
1152  w(4 ,iw) = ( dr + 5 )**9/( 2*256*dr6*( dr + 7 ) )
1153  w(5 ,iw) = -( dr + 7 )**9/( 256*dr8 )
1154  w(6 ,iw) = -( dr + 7 )**9/( 256*dr8 )
1155  w(7 ,iw) = -( dr + 7 )**9/( 256*dr8 )
1156  w(8 ,iw) = ( dr + 9 )**9/( 256*dr8*( dr + 9 ) )
1157  w(9 ,iw) = ( dr + 9 )**9/( 256*dr8*( dr + 9 ) )
1158  w(10,iw) = ( dr + 9 )**9/( 256*dr8*( dr + 9 ) )
1159  w(11,iw) = ( dr + 9 )**9/( 256*dr8*( dr + 9 ) )
1160  IF ( nd .GT. 2 ) THEN
1161  g(0,12) = 3/( dr + 9 )
1162  g(1,12) = 3/( dr + 9 )
1163  g(2,12) = 3/( dr + 9 )
1164  g(3,12) = 3/( dr + 9 )
1165  CALL smpcpy( 4, mxg, g(0,12), 1/( dr + 9 ) )
1166  pts(12) = ( ( dr + 1 )*dr*( dr - 1 )*( dr - 2 ) )/24
1167  w(12,iw) = w(8,iw)
1168  END IF
1169  END IF
1170 *
1171 * Compute constant weight values.
1172 *
1173  DO j = 1, rls
1174  w(1,j) = 1
1175  DO i = 2, wts
1176  w(1,j) = w(1,j) - pts(i)*w(i,j)
1177  END DO
1178  END DO
1179 *
1180 * Compute final weight values; null rule weights are computed as
1181 * differences between weights from highest degree and lower degree rules.
1182 *
1183  DO j = 2, rls
1184  DO i = 1, wts
1185  w(i,j) = w(i,j) - w(i,1)
1186  END DO
1187  END DO
1188 *
1189 ****END SMPRMS
1190 *
1191  END
1192 *
1193  SUBROUTINE smpcpy( START, END, PARAM, VALUE )
1194  DOUBLE PRECISION VALUE, PARAM(0:*)
1195  INTEGER START, END, I
1196  DO I = start, END
1197  param(i) = VALUE
1198  END DO
1199  END
1200 *
1201  SUBROUTINE smpsms( N, VERTEX, NF, F, G, SYMSMS, X, FUNVLS )
1202 *
1203 ****BEGIN PROLOGUE SMPSMS
1204 ****KEYWORDS fully symmetric sum
1205 ****PURPOSE To compute fully symmetric basic rule sums
1206 ****AUTHOR
1207 *
1208 * Alan Genz
1209 * Department of Mathematics
1210 * Washington State University
1211 * Pullman, WA 99164-3113, USA
1212 *
1213 ****LAST MODIFICATION 97-04
1214 ****DESCRIPTION SMPSMS computes a fully symmetric sum for a vector
1215 * of integrand values over a simplex. The sum is taken over
1216 * all permutations of the generators for the sum.
1217 *
1218 * ON ENTRY
1219 *
1220 * N Integer, number of variables.
1221 * VERTEX Real array of dimension (N,0:N)
1222 * The vertices of the simplex, one vertex per column.
1223 * NF Integer, number of components for the vector integrand.
1224 * F Subroutine for computing components of the integrand at X.
1225 * It must have parameters ( N, X, NF, FUNVLS );
1226 * Input parameters:
1227 * N Integer dimension of integral.
1228 * X Real array of length N, the evaluation point.
1229 * NF Integer number of components of the integrand.
1230 * Output parameter:
1231 * FUNVLS Real array of length NF, the integrand values at X.
1232 * G Real Array of dimension (0:N).
1233 * The generators for the fully symmetric sum.
1234 *
1235 * ON RETURN
1236 *
1237 * SYMSMS Real array of length NF, the values for the fully symmetric
1238 * sums for each component of the integrand.
1239 *
1240 ****ROUTINES CALLED: Integrand
1241 *
1242 ****END PROLOGUE SMPSMS
1243 *
1244 * Global variables.
1245 *
1246  INTEGER N, NF
1247  DOUBLE PRECISION VERTEX(N,0:N),G(0:N), SYMSMS(NF),FUNVLS(NF), X(N)
1248 *
1249 * Local variables.
1250 *
1251  INTEGER IX, LX, I, J, K, L
1252  DOUBLE PRECISION GL, GI
1253 *
1254 ****FIRST PROCESSING STATEMENT SymSum
1255 *
1256  DO i = 1, nf
1257  symsms(i) = 0
1258  END DO
1259 *
1260 * Sort generators if necessary
1261 *
1262  k = 0
1263  DO i = 1, n
1264  IF ( g(i) .GT. g(i-1) ) k = 1
1265  END DO
1266  IF ( k .GT. 0 ) THEN
1267  DO i = 1, n
1268  k = i - 1
1269  DO j = i, n
1270  IF ( g(j) .GT. g(k) ) k = j
1271  END DO
1272  IF ( k .GE. i ) THEN
1273  gi = g(i-1)
1274  g(i-1) = g(k)
1275  g(k) = gi
1276  END IF
1277  END DO
1278  END IF
1279 *
1280 * Compute integrand value for permutations of G
1281 *
1282  10 DO i = 1, n
1283  x(i) = vertex(i,0)*g(0)
1284  DO j = 1, n
1285  x(i) = x(i) + vertex(i,j)*g(j)
1286  END DO
1287  END DO
1288  CALL f( n, x, nf, funvls )
1289  DO j = 1, nf
1290  symsms(j) = symsms(j) + funvls(j)
1291  END DO
1292 *
1293 * Find next distinct permuation of G and loop back for value.
1294 * Permutations are generated in reverse lexicographic order.
1295 *
1296  DO i = 1, n
1297  IF ( g(i-1) .GT. g(i) ) THEN
1298  gi = g(i)
1299  ix = i - 1
1300  DO l = 0, i/2-1
1301  gl = g(l)
1302  g(l) = g(i-l-1)
1303  g(i-l-1) = gl
1304  IF ( gl .LE. gi ) ix = ix - 1
1305  IF ( g(l) .GT. gi ) lx = l
1306  END DO
1307  IF ( g(ix) .LE. gi ) ix = lx
1308  g(i) = g(ix)
1309  g(ix) = gi
1310  GO TO 10
1311  END IF
1312  END DO
1313 *
1314 ****END SMPSMS
1315 *
1316  END
1317 *
1318  SUBROUTINE smpstr( POINTR, SBRGNS, PONTRS, RGNERS )
1319 *
1320 ****BEGIN PROLOGUE SMPSTR
1321 ****AUTHOR
1322 *
1323 * Alan Genz
1324 * Department of Mathematics
1325 * Washington State University
1326 * Pullman, WA 99164-3113, USA
1327 *
1328 ****LAST MODIFICATION 2001-07
1329 ****PURPOSE SMPSTR maintains a heap for subregions.
1330 ****DESCRIPTION SMPSTR maintains a heap for subregions.
1331 * The subregions are ordered according to the size of the
1332 * greatest error estimates of each subregion (RGNERS).
1333 *
1334 * PARAMETERS
1335 *
1336 * POINTR Integer.
1337 * The index for the subregion to be inserted in the heap.
1338 * SBRGNS Integer.
1339 * Number of subregions in the heap.
1340 * PONTRS Real array of dimension SBRGNS.
1341 * Used to store the indices for the greatest estimated errors
1342 * for each subregion.
1343 * RGNERS Real array of dimension SBRGNS.
1344 * Used to store the greatest estimated errors for each
1345 * subregion.
1346 *
1347 ****ROUTINES CALLED NONE
1348 ****END PROLOGUE SMPSTR
1349 *
1350 * Global variables.
1351 *
1352  INTEGER POINTR, SBRGNS
1353  DOUBLE PRECISION PONTRS(*), RGNERS(*)
1354 *
1355 * Local variables.
1356 *
1357 * RGNERR Intermediate storage for the greatest error of a subregion.
1358 * SUBRGN Position of child/parent subregion in the heap.
1359 * SUBTMP Position of parent/child subregion in the heap.
1360 *
1361  INTEGER SUBRGN, SUBTMP, PT, PTP
1362  DOUBLE PRECISION RGNERR
1363 *
1364 ****FIRST PROCESSING STATEMENT SMPSTR
1365 *
1366  RGNERR = rgners(pointr)
1367  IF ( pointr .EQ. pontrs(1) ) THEN
1368 *
1369 * Move the new subregion inserted at the top of the heap
1370 * to its correct position in the heap.
1371 *
1372  subrgn = 1
1373  10 subtmp = 2*subrgn
1374  IF ( subtmp .LE. sbrgns ) THEN
1375  IF ( subtmp .NE. sbrgns ) THEN
1376 *
1377 * Find maximum of left and right child.
1378 *
1379  pt = pontrs(subtmp)
1380  ptp = pontrs(subtmp+1)
1381  IF ( rgners(pt) .LT. rgners(ptp) ) subtmp = subtmp + 1
1382  ENDIF
1383 *
1384 * Compare maximum child with parent.
1385 * If parent is maximum, then done.
1386 *
1387  pt = pontrs(subtmp)
1388  IF ( rgnerr .LT. rgners(pt) ) THEN
1389 *
1390 * Move the pointer at position subtmp up the heap.
1391 *
1392  pontrs(subrgn) = pt
1393  subrgn = subtmp
1394  GO TO 10
1395  ENDIF
1396  ENDIF
1397  ELSE
1398 *
1399 * Insert new subregion in the heap.
1400 *
1401  subrgn = sbrgns
1402  20 subtmp = subrgn/2
1403  IF ( subtmp .GE. 1 ) THEN
1404 *
1405 * Compare child with parent. If parent is maximum, then done.
1406 *
1407  pt = pontrs(subtmp)
1408  IF ( rgnerr .GT. rgners(pt) ) THEN
1409 *
1410 * Move the pointer at position subtmp down the heap.
1411 *
1412  pontrs(subrgn) = pt
1413  subrgn = subtmp
1414  GO TO 20
1415  ENDIF
1416  ENDIF
1417  ENDIF
1418  pontrs(subrgn) = pointr
1419 *
1420 ****END SMPSTR
1421 *
1422  END
1423 *
1424  SUBROUTINE smpdfs( ND, NF, FUNSUB, TOP, SBRGNS, VERTCS, VOLUMS,
1425  * X, H, CENTER, WORK, FRTHDF, NEWSBS )
1426 *
1427 ****BEGIN PROLOGUE SMPDFS
1428 ****PURPOSE To compute new subregions
1429 ****AUTHOR
1430 *
1431 * Alan Genz
1432 * Department of Mathematics
1433 * Washington State University
1434 * Pullman, WA 99164-3113, USA
1435 *
1436 ****LAST MODIFICATION 2001-07
1437 ****DESCRIPTION SMPDFS computes fourth differences along each edge
1438 * direction. It uses these differences to determine a
1439 * subdivision of the orginal subregion into either three or
1440 * four new subregions.
1441 *
1442 * ON ENTRY
1443 *
1444 * ND Integer, number of variables.
1445 * NF Integer, number of components for the vector integrand.
1446 * FUNSUB Externally declared subroutine.
1447 * For computing the components of the integrand at a point X.
1448 * It must have parameters (ND, X, NF, FUNVLS).
1449 * Input Parameters:
1450 * X Real array of dimension ND, the evaluation point.
1451 * ND Integer, number of variables for the integrand.
1452 * NF Integer, number of components for the vector integrand.
1453 * Output Parameters:
1454 * FUNVLS Real array of dimension NF.
1455 * The components of the integrand at the point X.
1456 * TOP Integer, location in VERTCS array for original subregion.
1457 * SBRGNS Integer, number of subregions in VERTCS BEFORE subdivision.
1458 * VERTCS Real array of dimension (ND,0:ND,*), vertices of orginal
1459 * subregion must be in VERTCS(1:ND,0:ND,TOP).
1460 * VOLUMS Real array of dimension (*) of volumes for subregions.
1461 * X Real work array of dimension ND.
1462 * H Real work array of dimension ND.
1463 * CENTER Real work array of dimension (0:ND).
1464 * WORK Real work array of dimension 5*NF.
1465 * FRTHDF Real work array of dimension (0:ND-1,ND).
1466 *
1467 * ON RETURN
1468 *
1469 * NEWSBS Integer, number of new subregions (3 or 4).
1470 * FUNCLS Integer, number of FUNSUB calls used by SMPDFS.
1471 * VERTCS Real array of dimension (ND,0:ND,*).
1472 * The vertices of the of new subegions will be at locations
1473 * TOP, SBRGNS+1, ..., SBRGNS+NEWSBS-1.
1474 * VOLUMS Real Array of dimension (*).
1475 * VOLUMS has been updated for new subregions.
1476 *
1477 ****ROUTINES CALLED: FUNSUB
1478 *
1479 ****END PROLOGUE SMPDFS
1480 *
1481  EXTERNAL funsub
1482  INTEGER ND, NF, TOP, SBRGNS, NEWSBS
1483  DOUBLE PRECISION VERTCS(ND,0:ND,*), VOLUMS(*), WORK(NF,*)
1484  DOUBLE PRECISION X(ND), H(ND), CENTER(ND), FRTHDF(0:ND-1,ND)
1485  DOUBLE PRECISION DIFFER, DIFMAX, DIFMID, DIFNXT, EWIDTH, EDGMAX
1486  DOUBLE PRECISION CUTTF, CUTTB, DIFIL, DIFLJ, DFSMAX, VTI, VTJ, VTL
1487  PARAMETER ( CUTTF = 2, cuttb = 8 )
1488  INTEGER I, J, K, L, IE, JE, IS, JS, LS, IT, JT, LT
1489  DOUBLE PRECISION SMPVOL
1490 *
1491 ****FIRST PROCESSING STATEMENT SMPDFS
1492 *
1493 *
1494 * Compute the differences.
1495 *
1496  is = 0
1497  js = 1
1498  difmax = 0
1499  edgmax = 0
1500  DO k = 1, nd
1501  center(k) = vertcs(k,0,top)
1502  DO l = 1, nd
1503  center(k) = center(k) + vertcs(k,l,top)
1504  END DO
1505  center(k) = center(k)/( nd + 1 )
1506  END DO
1507  CALL funsub(nd, center, nf, work(1,3))
1508  DO i = 0, nd-1
1509  DO j = i+1, nd
1510  ewidth = 0
1511  DO k = 1, nd
1512  h(k) = 2*( vertcs(k,i,top)-vertcs(k,j,top) )/( 5*(nd+1) )
1513  ewidth = ewidth + abs( h(k) )
1514  x(k) = center(k) - 3*h(k)
1515  END DO
1516  DO l = 1, 5
1517  DO k = 1, nd
1518  x(k) = x(k) + h(k)
1519  END DO
1520  IF ( l. ne. 3 ) CALL funsub(nd, x, nf, work(1,l))
1521  END DO
1522  IF ( ewidth .GE. edgmax ) THEN
1523  ie = i
1524  je = j
1525  edgmax = ewidth
1526  END IF
1527  differ = 0
1528  difmid = 0
1529  DO k = 1, nf
1530  difmid = difmid + abs( work(k,3) )
1531  differ = differ + abs( work(k,1) + work(k,5)+ 6*work(k,3)
1532  & - 4*( work(k,2) + work(k,4) ) )
1533  END DO
1534  IF ( difmid + differ/8 .EQ. difmid ) differ = 0
1535  differ = differ*ewidth
1536  frthdf(i,j) = differ
1537  IF ( differ .GE. difmax ) THEN
1538  it = is
1539  jt = js
1540  difnxt = difmax
1541  is = i
1542  js = j
1543  difmax = differ
1544  ELSE IF ( differ .GE. difnxt ) THEN
1545  it = i
1546  jt = j
1547  difnxt = differ
1548  END IF
1549  END DO
1550  END DO
1551 *
1552 * Determine whether to compute three or four new subregions.
1553 *
1554  IF ( difnxt .GT. difmax/cuttf ) THEN
1555  newsbs = 4
1556  ELSE
1557  newsbs = 3
1558  IF ( difmax .EQ. 0 ) THEN
1559  is = ie
1560  js = je
1561  ELSE
1562  dfsmax = 0
1563  DO l = 0, nd
1564  IF ( l .NE. is .AND. l .NE. js ) THEN
1565  it = min( l, is, js )
1566  jt = max( l, is, js )
1567  lt = is + js + l - it - jt
1568  differ = frthdf(it,lt) + frthdf(lt,jt)
1569  IF ( differ .GE. dfsmax ) THEN
1570  dfsmax = differ
1571  ls = l
1572  END IF
1573  END IF
1574  END DO
1575  difil = frthdf( min(is,ls), max(is,ls) )
1576  diflj = frthdf( min(js,ls), max(js,ls) )
1577  difnxt = difil + diflj - min( difil,diflj )
1578  IF ( difmax/cuttb .LT. difnxt .AND. difil .GT. diflj ) THEN
1579  it = is
1580  is = js
1581  js = it
1582  END IF
1583  END IF
1584  END IF
1585 *
1586 * Copy vertices and volume for TOP to new subregions
1587 *
1588  volums(top) = volums(top)/newsbs
1589  DO l = sbrgns + 1, sbrgns + newsbs - 1
1590  volums(l) = volums(top)
1591  DO j = 0, nd
1592  DO k = 1, nd
1593  vertcs(k,j,l) = vertcs(k,j,top)
1594  END DO
1595  END DO
1596  END DO
1597  DO k = 1, nd
1598  vti = vertcs(k,is,top)
1599  vtj = vertcs(k,js,top)
1600  IF ( newsbs .EQ. 4 ) THEN
1601 *
1602 * Compute four new subregions.
1603 *
1604  vertcs(k,js,top) = ( vti + vtj )/2
1605  vertcs(k,is,sbrgns+1) = vti
1606  vertcs(k,js,sbrgns+1) = ( vti + vtj )/2
1607  vertcs(k,is,sbrgns+2) = ( vti + vtj )/2
1608  vertcs(k,js,sbrgns+2) = vtj
1609  vertcs(k,is,sbrgns+3) = ( vti + vtj )/2
1610  vertcs(k,js,sbrgns+3) = vtj
1611  vti = vertcs(k,it,top)
1612  vtj = vertcs(k,jt,top)
1613  vertcs(k,jt,top) = ( vti + vtj )/2
1614  vertcs(k,it,sbrgns+1) = ( vti + vtj )/2
1615  vertcs(k,jt,sbrgns+1) = vtj
1616  vti = vertcs(k,it,sbrgns+2)
1617  vtj = vertcs(k,jt,sbrgns+2)
1618  vertcs(k,jt,sbrgns+2) = ( vti + vtj )/2
1619  vertcs(k,it,sbrgns+3) = ( vti + vtj )/2
1620  vertcs(k,jt,sbrgns+3) = vtj
1621  ELSE
1622 *
1623 * Compute three new subregions.
1624 *
1625  vertcs(k,js,top) = ( 2*vti + vtj )/3
1626  vertcs(k,is,sbrgns+1) = ( 2*vti + vtj )/3
1627  IF ( difmax/cuttf .LT. difnxt ) THEN
1628  vertcs(k,js,sbrgns+1) = vtj
1629  vertcs(k,is,sbrgns+2) = ( 2*vti + vtj )/3
1630  vertcs(k,js,sbrgns+2) = vtj
1631  vtj = vertcs(k,js,sbrgns+1)
1632  vtl = vertcs(k,ls,sbrgns+1)
1633  vertcs(k,ls,sbrgns+1) = ( vtj + vtl )/2
1634  vertcs(k,js,sbrgns+2) = ( vtj + vtl )/2
1635  vertcs(k,ls,sbrgns+2) = vtl
1636  ELSE
1637  vertcs(k,js,sbrgns+1) = ( vti + 2*vtj )/3
1638  vertcs(k,is,sbrgns+2) = ( vti + 2*vtj )/3
1639  vertcs(k,js,sbrgns+2) = vtj
1640  END IF
1641  END IF
1642  END DO
1643 *
1644 ****END SMPDFS
1645 *
1646  END
1647 
1648 
vertex
A vertex.
Definition: vertex.hpp:21