* ************************************************************************* subroutine fdelgf( iel, elptr, elvar, xst, x, + gptr, fuval, nefval, inform ) * ************************************************************************* * Purpose : * --------- * This routine estimates the IELth element gradient vector * by forward differences in the element function values. * This IELth element gradient has an elemental representation. * Every time information on the element function value is needed, * a return is made to the main program. This routine is then * re-entered with the information available and the calculation * proceeds. * Parameters : * ------------ * iel ( int ) * input : the indice of the element gradient vector * that will be estimated. * output : unmodified. * elptr ( int ) * input : array whose kth value is the position of * the first variable of element k, in the * list ELVAR. * output : unmodified. * elvar ( int ) * input : array containing the indices of the varaibles * in the first element, followed by those in the * second element, etc. * output : unmodified. * xst ( int ) * input : vector containing the status of the variables. * output : unmodified. * x ( dble ) * input : the current iterate vector. * output : unmodified. * gptr ( int ) * input : array whose ith value is the position of the * first component of the ith element gradient * in FUVAL. * output : unmodified. * fuval ( dble ) * input : array used to store the function and derivative * values of the element functions. * output : the components of indices GPTR(IEL) up to * GPTR(IEL+1)-1 will be set to the estimate * of the element gradient vector by forward * differences in the element function values. * nefval ( int ) * input : the number of element function evaluations. * output : this number is increased by one each time * one element function needs to be evaluated. * inform ( int ) * input : must be set to zero the first time this routine * is called. * If it is not equal to zero, the routine is re-entered * with an element function value that was needed * to proceed the calculation. * output : If it is not equal to zero, an element function * value is needed to proceed the calculation. A return * is then made to the main program. * Otherwise, no further information is needed. * Routines used : * --------------- * sqrt, max, abs, sign. * gxfix. * Programming : * ------------- * D. Tuyttens * ======================================================================== * Routine parameters integer iel, elptr(*), elvar(*), xst(*), + gptr(*), nefval, inform double precision x(*), fuval(*) * Internal variables integer i, k, ik, debg, istr, iend logical gxfix double precision sqrem, step, tempx, tempf * Save variables save i, k, ik, debg, istr, iend save sqrem, step, tempx, tempf * Common specifications double precision epsmch, huge, tiny, tol common / prbmch / epsmch, huge, tiny, tol double precision zero, one, two, three, half, tenm1, tenm2, tenm4 common / prbcst / zero, one, two, three, half, tenm1, tenm2, tenm4 * * Reverse communication test. * if( inform.ne.0 ) goto 100 * * Estimate of the IELth element gradient vector by * forward differences in the element function values. * sqrem = sqrt(epsmch) debg = gptr(iel)-1 tempf = fuval(iel) ik = 0 istr = elptr(iel) iend = elptr(iel+1)-1 if ( iend .ge. istr ) then i = istr 10 continue k = elvar(i) ik = ik + 1 if( gxfix(xst(k)) ) then fuval(debg+ik) = zero else tempx = x(k) step = sqrem * sign( max(abs(tempx),one) , tempx ) x(k) = tempx + step nefval = nefval + 1 inform = 6 return * * Element function evaluation. * 100 continue inform = 0 x(k) = tempx fuval(debg+ik) = ( fuval(iel) - tempf ) / step endif i = i + 1 if ( i .le. iend ) go to 10 endif fuval(iel) = tempf * return end