Lib-DVM interface description (contents) Part 1
(1-5)
Part 2
(6-7)
Part 3
(8-11)
Part 4
(12-13)
Part 5
(14-15)
Part 6
(16-18)
Part 7
(19)
created: february, 2001 - last edited 03.05.01 -

19 Examples of programs using Run-Time System functions

19.1 Solution of Laplas equation by Jacobi method

The usage of main functions of Run-Time System is shown on the example of numerical solution of two-dimensional Laplas equation in specified rectangular area with specified boundary conditions (Dirixle task)

2U
¾ ¾ +
x2
2U
¾ ¾ = 0 (1)
y2

Three-point approximation of partial derivatives of the equation (1) on squared grid using five-point pattern

(i-1,j) , (i+1,j) , (i,j) , (i,j-1) , (i,j+1)

gives differential equation ("cross" scheme), reduced to the equation

ui,j = (ui+1,j + ui-1,j + ui,j+1 + ui,j-1)/4 (i,j=1, ... , k-2) (2) ,

where:

k - a number of grid points of each dimension;
ui,j - a value of grid function in a node (i,j).

The values of the grid function

u0,j , uk-1,j , ui,0 , ui,k-1 (i,j=0, ... , k-1)

are determined by boundary conditions.

In the programs below (Fortran and C) the solution of the equation (1) is seached by formulas (2) using method of sequential approximations. The iterational proccess is completed if maximal diviation of the grid function ui,j on all knots for to sequential iterations will become less, then a value of variable maxeps. The initial values of the grid functoin are ui,j = 1+i+j (i,j = 1, … , k-2). The boundary conditions are zero: u0,j=0, uk-1,j=0, ui,0=0, ui,k-1=0 (i,j = 0, ... , k-1). The solution is searched on square [0:7, 0:7].

PROGRAM IN FORTRAM LANGUAGE

    program cross
      integer linit, lexit, getam, getps, crtamv, distr, crtda, align,
     +        crtpl, dvmadr, mappl, endpl, exfrst, imlast, dopl,
     +        tstio, crtrg, crtred, insred, strtrd, waitrd,
     +        crtshg, inssh, strtsh, waitsh, delrg, delshg
      real bptr(1)
      integer dvm
      integer amref, psref, mvref, plref, rgref, redref, shgref
      integer amdim(2), disaxs(2), dispar(2)
      integer shwdth(2), axis(2), coeff(2), const(2)
      integer lvaddr(2), lvtype(2), iiniti(2), ilasti(2), istep(2)
      integer oiniti(2), olasti(2), ostep(2)
C     A number of grid points by each dimension and
C     maximal number of iterations
      parameter (k = 8, itmax = 20)
      real eps, maxeps
C     Array header with previous values of the grid function
      integer ahdr(3)
C     Array header with next values of the grid function
      integer bhdr(3)
      maxeps = 0.5e-7
C     Run-Time System initialisation
      dvm = linit (0)
C     Creating abstract machine representation
C     and mapping it onto processor subsystem
      amref = getam ()
      psref = getps (amref)
      amdim(1) = k
      amdim(2) = k
      mvref = crtamv (amref, 2, amdim,0)
      disaxs(1) = 1
      disaxs(2) = 2
      dispar(1) = 0
      dispar(2) = 0
      dvm = distr (mvref, psref, 2, disaxs, dispar)
C     Creating and mapping arrays
C     with the values of the grid function
      shwdth(1) = 1
      shwdth(2) = 1
      dvm = crtda (ahdr, 1, bptr, 2, 4, amdim, 0, 0, shwdth, shwdth)
      dvm = crtda (bhdr, 1, bptr, 2, 4, amdim, 0, 0, shwdth, shwdth)
      axis(1) = 1
      axis(2) = 2
      coeff(1) = 1
      coeff(2) = 1
      const(1) = 0
      const(2) = 0
      dvm = align (ahdr, mvref, axis, coeff, const)
      dvm = align (bhdr, mvref, axis, coeff, const)
C     Parallel loop of initializing arrays
C     with values of the grid function
C     (parallel loop with the base array ahdr)
      plref = crtpl (2)
      lvaddr(1) = dvmadr (j)
      lvaddr(2) = dvmadr (i)
      lvtype(1) = 1
      lvtype(2) = 1 
      iiniti(1) = 0
      iiniti(2) = 0
      ilasti(1) = k - 1
      ilasti(2) = k - 1
      istep(1) = 1
      istep(2) = 1
      dvm = mappl (plref, ahdr, axis, coeff, const, lvaddr, lvtype,
     +             iiniti, ilasti, istep, oiniti, olasti, ostep)
4 if (dopl (plref) .eq. 0) goto 5
      do 1 j = oiniti(1), olasti(1), ostep(1)
         do 1 i = oiniti(2), olasti(2), ostep(2)
            bptr( ahdr(3) + 1 + i + ahdr(2) * j ) = 0.
            bptr( bhdr(3) + 1 + i + bhdr(2) * j ) = 1. + i + j
1     continue
      goto 4
5     dvm = endpl (plref)
C     Creating reduction variable and reduction group
C     to calculate maximal deviation of the grid function for
C     two sequential iterations
      redref = crtred ( 3, eps, 3, 1, 0, 0, 0)
      rgref = crtrg (0, 0)
      dvm = insred (rgref, redref, 0)
C     Creating shadow edge group for renewing shadow edges
      shgref = crtshg (0)
      dvm = inssh (shgref, ahdr, shwdth, shwdth, 0)
C     MAIN ITERATION LOOP
      do 2 it = 1,itmax
C     Parallel loop to calculate maximal deviation
C     of the grid function in variable eps
C     (parallel loop with base array ahdr)
      plref = crtpl (2)
      iiniti(1) = 1
      iiniti(2) = 1
      ilasti(1) = k - 2
      ilasti(2) = k - 2
      dvm = mappl (plref, ahdr, axis, coeff, const, lvaddr, lvtype,
    +              iiniti, ilasti, istep, oiniti, olasti, ostep)
      eps = 0.
6     if (dopl (plref) .eq. 0) goto 8
      do 7 j = oiniti(1), olasti(1), ostep(1)
         do 7 i = oiniti(2), olasti(2), ostep(2)
            eps = max(eps, abs(bptr( bhdr(3)+1+i+bhdr(2)*j )
     +                       - bptr( ahdr(3)+1+i+ahdr(2)*j )))
            bptr( ahdr(3) + 1 + i + ahdr(2) * j ) = 
     +      bptr( bhdr(3) + 1 + i + bhdr(2) * j )
7     continue
      goto 6
8     dvm = endpl (plref)
C     Calculation of reductional maximum
      dvm = strtrd (rgref)
      dvm = waitrd (rgref)
C     Shadow edge exchange
      dvm = strtsh (shgref)
      dvm = waitsh (shgref)
C     Parallel loop to calculate new values
C     of the grid function in array bhdr
C     (parallel loop with base array bhdr)
      plref = crtpl (2)
      dvm = mappl ( plref, bhdr, axis, coeff, const, lvaddr, lvtype,
     +              iiniti, ilasti, istep, oiniti, olasti, ostep)
9     if (dopl (plref) .eq. 0) goto 11
      do 10 j = oiniti(1), olasti(1), ostep(1)
         do 10 i = oiniti(2), olasti(2), ostep(2)
              bptr(bhdr(3) + 1 + i + bhdr(2) * j) =
     +        ( bptr(ahdr(3) + 1 + (i - 1) + ahdr(2) * j) +
     +          bptr(ahdr(3) + 1 + i + ahdr(2) * (j - 1)) +
     +          bptr(ahdr(3) + 1 + (i + 1) + ahdr(2) * j) +
     +          bptr(ahdr(3) + 1 + i + ahdr(2) * (j + 1)) ) / 4.
10    continue
      goto 9
11    dvm = endpl (plref)
C     Printing current deviation of the grid function
      if (tstio ()) print *,'IT = ',it,' EPS = ',eps
C     Exit, if specified precision is achieved
      if (eps .lt. maxeps) goto 3
2     continue
C     END OF MAIN ITERATION LOOP
C     End of dealing with Run-Time System
3     dvm = lexit (0)
      end

PROGRAM IN C LANGUAGE

#include "dvmlib.h"
int  main(int argc, char *argv[])
{ 
   long  InitPar = 0, UserRes = 0, ExtHdr = 0;
   long  StaticMV = 0, StaticDA = 0, ReDistrDA = 0, StaticRV = 0,
         DelRed = 0, StaticRG = 0, StaticShG = 0, FullShd = 0;
   long  i, j, it, Rank2 = 2, Long0 = 0, TypeSize = sizeof(float);
   long  RedFunc = rf_MAX, RVType = rt_FLOAT, RedArrayLength = 1;
         AMRef amref;
   PSRef           psref;
   AMViewRef       mvref;
   LoopRef         plref;
   RedGroupRef     rgref;
   RedRef          redref;
   ShadowGroupRef  shgref;
   AddrType        lvaddr[2];
   long            lvtype[2] = {0,0};
   long            amdim[2], disaxs[2], dispar[2];
   long            shwdth[2], axis[2], coeff[2], cnst[2];
   long            iiniti[2], ilasti[2], istep[2];
   long            oiniti[2], olasti[2], ostep[2];
   DVMFILE        *OutFile;
   float  eps;                 /* current precision of calculation */
    long  k = 8;               /* a number of grid points by every dimension */
   float maxeps = 0.5e-2f; /* requied precision of calculation */
   long itmax = 200;       /* maximal number of iterations */
   long Ahdr[3];           /* header of array with previous 
                              values of the grid function */
   long Bhdr[3];           /* header of array with next
                              values of the grid function */
   rtl_init(InitPar, argc, argv); /* Run-Time System initialization */
   /* Creating abstract machine representation
      and its mapping on processor subsystem */
   amref = getam_ ();
   psref = getps_ (&amref);
   amdim[0] = k;
   amdim[1] = k;
   mvref = crtamv_ (&amref, &Rank2, amdim, &StaticMV);
   disaxs[0] = 1;
   disaxs[1] = 2;
   dispar[0] = 0;
   dispar[1] = 0;
   distr_ (&mvref, &psref, &Rank2, disaxs, dispar);
   /* Creation and distribution of arrays
      with values of the grid function */
   shwdth[0] = 1;
   shwdth[1] = 1;
   /* Array with previous values of the grig function */
   crtda_ (Ahdr, &ExtHdr, NULL, &Rank2, &TypeSize, amdim,
           &StaticDA, &ReDistrDA,shwdth, shwdth);
   /* Array with next values of grid function */
   crtda_ (Bhdr, &ExtHdr, NULL, &Rank2, &TypeSize, amdim,
           &StaticDA, &ReDistrDA, shwdth, shwdth);
   axis[0] = 1;
   axis[1] = 2;
   coeff[0] = 1;
   coeff[1] = 1;
   cnst[0] = 0;
   cnst[1] = 0;
   align_ (Ahdr,&mvref,axis,coeff,cnst);
   align_ (Bhdr,&mvref,axis,coeff,cnst);
   /* Parallel loop of initializing arrays
      with values of the grid function
      (parallel loop with base array Ahdr) */
   plref = crtpl_ (&Rank2);
   lvaddr[0] = (AddrType)&j;
   lvaddr[1] = (AddrType)&i;
   iiniti[0] = 0;
   iiniti[1] = 0;
   ilasti[0] = k - 1;
   ilasti[1] = k - 1;
   istep[0] = 1;
   istep[1] = 1;
   mappl_ (&plref, (PatternRef *)Ahdr, axis, coeff, cnst, lvaddr,
           lvtype, iiniti, ilasti, istep, oiniti, olasti, ostep);
   while (dopl_ (&plref))
      for ( j = oiniti[0]; j <= olasti[0]; j += ostep[0] )
         for ( i = oiniti[1]; i <= olasti[1]; i += ostep[1] )
         {  DAElm2(Ahdr, float, j, i) = 0.;
            DAElm2(Bhdr, float, j, i) = 1.f + i + j;
         }
   endpl_ (&plref);
   /* Creating reduction variable and reduction group
      to calculate maximal deviation of the grid function for
      two sequential iterations */
   redref = crtred_ (&RedFunc, &eps, &RVType, &RedArrayLength,
                     NULL, &Long0, &StaticRV);
   rgref = crtrg_ (&StaticRG, &DelRed);
   insred_ (&rgref, &redref, NULL);
   /* Creating shadow edge group for renewing shadow edges */
   shgref = crtshg_ (&StaticShG);
   inssh_ (&shgref, Ahdr, shwdth, shwdth, &FullShd);
   /* MAIN ITERATION LOOP */
   for (it = 1; it <= itmax; it++)
   {
       /* Parallel loop to calculate maximal deviation
          of the grid function in variable eps
          (parallel loop with base array Ahdr) */
       plref = crtpl_ (&Rank2);
       iiniti[0] = 1;
       iiniti[1] = 1;
       ilasti[0] = k - 2;
       ilasti[1] = k - 2;
       mappl_ (&plref, (PatternRef *)Ahdr, axis, coeff, cnst, lvaddr,
               lvtype, iiniti, ilasti, istep, oiniti, olasti, ostep);
       eps = 0.;
       while (dopl_ (&plref))
          for ( j = oiniti[0]; j <= olasti[0]; j += ostep[0] )
             for ( i = oiniti[1]; i <= olasti[1]; i += ostep[1] )
             { eps = max(eps, dvm_abs( DAElm2(Bhdr,float,j,i) -
                         DAElm2(Ahdr,float,j,i) ));
               DAElm2(Ahdr, float, j, i) =
               DAElm2(Bhdr, float, j, i);
             } 
       endpl_ (&plref);
       /* Calculation of reductional maximum */
       strtrd_ (&rgref);
       waitrd_ (&rgref);
       /* Shadow edges exchange */
       strtsh_ (&shgref);
       waitsh_ (&shgref);
       /* Parallel loop to calculate new values
          of the grid function in array Bhdr
          (parallel loop with base array Bhdr) */
       plref = crtpl_ (&Rank2);
       mappl_ (&plref, (PatternRef *)Bhdr, axis, coeff, cnst, lvaddr,
               lvtype, iiniti, ilasti, istep, oiniti, olasti, ostep);
       while (dopl_ (&plref))
          for ( j = oiniti[0]; j <= olasti[0]; j += ostep[0] )
             for ( i = oiniti[1]; i <= olasti[1]; i += ostep[1] )
                 DAElm2(Bhdr, float, j, i) =
                 ( DAElm2(Ahdr, float, j, i-1) +
                   DAElm2(Ahdr, float, j-1, i) +
                   DAElm2(Ahdr, float, j, i+1) +
                   DAElm2(Ahdr, float, j+1, i) ) / 4.f;
       endpl_ (&plref);
       /* Printing current deviation of grid function */
       dvm_printf("IT = %ld EPS = %e\n", it, eps);
       if (eps < maxeps)
           break; /* exit, if specified precision is achieved */  
   }
   /* END OF MAIN ITERATION LOOP */
   /* Writing solution in the file
      and end of work with Run-Time System */
   OutFile = dvm_fopen("cross.dat", "w+b");
   if (OutFile)
      { dvm_dfwrite(Bhdr, 0, OutFile);
        dvm_fclose(OutFile);
      }
   lexit_ (&UserRes);
   return (int)UserRes;
}

19.2 Parallel loop with regular data dependence between iterations

The (C and Fortran) programs, considered below, are the examples of implementation of parallel loop with data (flow) dependence and anti-dependence between iterations using Run-Time System functions. The parallel loop is two-dimensional loop of the form

for( j = 1; j <= k-2; j++ )  
  for( i = 1; i <= k-2; i++ )  
    A[j,i] = ( A[j- FDL1, i] + /* flow dependence on dimension 1 */
    A[j+ ADL1, i] + /* anti-dependence on dimension 1 */
    A[j, i- FDL2] + /* flow dependence on dimension 2 */
    A[j, i+ ADL2] /* anti-dependence on dimension 2 */
    ) / 4;  

Where:

k size of A array dimensions;
FDL1 dependence length for dimension 1;
ADL1 anti-dependence length for dimension 1;
FDL2 dependence length for dimension 2;
ADL2 anti-dependence length for dimension 2.

PROGRAM IN FORTRAM LANGUAGE

     program across
      integer linit, lexit,getam, getps, crtamv, distr, crtda, align,
     +        crtpl, dvmadr, mappl, endpl, dopl,
     +        crtshg, inssh, insshd, strtsh, recvsh, sendsh, waitsh
      real    bptr(1)
      integer dvm
      integer amref, psref, mvref, plref, shgref
      integer amdim(2), disaxs(2), dispar(2)
      integer lshwd(2), hshwd(2), shsign(2), axis(2), coeff(2),
     +        const(2)
      integer lvaddr(2), lvtype(2), iiniti(2), ilasti(2), istep(2)
      integer oiniti(2), olasti(2), ostep(2)
      integer FDL1, FDL2, ADL1, ADL2, finssh
C     Size of each array dimension
C     with regular data dependence
      parameter (k = 200)
C     Header of array with regular data dependence
      integer ahdr(6)
      ahdr(5) = 1
      ahdr(6) = 1
C     Length of data dependence for dimension 1
      FDL1 = 4
C     Length of data dependence for dimension 2
      FDL2 = 3
C     Length of data anti-dependence for dimension 1
      ADL1 = 2
C     Length of data anti-dependence for dimension 2
      ADL2 = 1
C     Flag to use inssh_ function to registrate array
C     in shadow edge group;
C     If finssh=0, insshd_ function is used
      finssh = 0
C     Run-Time System initialization
      dvm = linit (1)
C     Creating abstract machine representation and
C     its mapping on processor subsystem
      amref = getam ()
      psref = getps (amref)
      amdim(1) = k
      amdim(2) = k
      mvref = crtamv (amref, 2, amdim,0)
      disaxs(1) = 1
      disaxs(2) = 2
      dispar(1) = 0
      dispar(2) = 0
      dvm = distr (mvref, psref, 2, disaxs, dispar)
C     Creation and distribution of array
C     with regular data dependence
C     Width of low shadow edge for dimension 1 is
C     length of data dependence for dimension 1
      lshwd(1) = FDL1
C     Width of low shadow edge for dimension 2 is
C     length of data dependence for dimension 2
      lshwd(2) = FDL2
C     Width of high shadow edge for dimension 1 is
C     length of data anti-dependence for dimension 1
      hshwd(1) = ADL1
C     Width of high shadow edge for dimension 2 is
C     length of data anti-dependence for dimension 2
      hshwd(2) = ADL2
      dvm = crtda (ahdr, 1, bptr, 2, 4, amdim, 0, 0, lshwd, hshwd)
      axis(1) = 1
      axis(2) = 2
      coeff(1) = 1
      coeff(2) = 1
      const(1) = 0
      const(2) = 0
      dvm = align (ahdr, mvref, axis, coeff, const)
C     Parallel loop of array initialization
C     with regular data dependence
      plref = crtpl (2)
      lvaddr(1) = dvmadr (j)
      lvaddr(2) = dvmadr (i)
      lvtype(1) = 1
      lvtype(2) = 1
      iiniti(1) = 0
      iiniti(2) = 0
      ilasti(1) = k - 1
      ilasti(2) = k - 1
      istep(1) = 1
      istep(2) = 1
      dvm = mappl (plref, ahdr, axis, coeff, const, lvaddr, lvtype
     +             iiniti, ilasti, istep, oiniti, olasti, ostep)
2     if (dopl (plref) .eq. 0) goto 3
      do 1 j = oiniti(1), olasti(1), ostep(1)
         do 1 i = oiniti(2), olasti(2), ostep(2)
            bptr( ahdr(3) + 1 + i + ahdr(2) * j ) = 1.
1     continue
      goto 2
3     dvm = endpl (plref)
C     Creation of shadow edge group to supply
C     data anti-dependence
      shgref = crtshg (0)
C     Zero width of low shadow edge for dimension 1
      lshwd(1) = 0
C     Zero width of low shadow edge for dimension 2
      lshwd(2) = 0
C     Width of high shadow edge for dimension 1 is
C     length of data anti-dependence for dimension 1
      hshwd(1) = ADL1
C     Width of high shadow edge for dimension 2 is
C     length of data anti-dependence for dimension 2
      hshwd(2) = ADL2
      dvm = inssh (shgref, ahdr, lshwd, hshwd, 0)
C     Initialization of shadow edges to support
C     data anti-dependence
      dvm = strtsh (shgref)
      dvm = waitsh (shgref)
C     Creation of shadow edge group to support
C     data dependence
      shgref = crtshg (0)
C     Width of low shadow edge for dimension 1 is
C     length of data dependence for dimension 1
      lshwd(1) = FDL1
C     Width of low shadow edge for dimension 2 is
C     length of data dependence for dimension 2
      lshwd(2) = FDL2
C     Zero width of high shadow edge for dimension 1
      hshwd(1) = 0
C     Zero width of high shadow edge for dimension 2
      hshwd(2) = 0 
      if (finssh .eq. 0) goto 4
C     inssh_ is used to registrate array in shadow edge group
      dvm = inssh (shgref, ahdr, lshwd, hshwd, 0)
      goto 5
C     insshd_ is used to registrate array in shadow edge group
C     Only low edge is used for dimensions 1 and 2
4     shsign(1) = 3
      shsign(2) = 3
      dvm = insshd (shgref, ahdr, lshwd, hshwd, 1, shsign)
      goto 5
C     Receiving imported elements
C     to support data dependence
5     dvm = recvsh (shgref)
      dvm = waitsh (shgref)
C     PARALLEL LOOP WITH DATA DEPENDENCE AND
C     ANTI-DEPENDENCE BETWEEN ITERATIONS
      plref = crtpl (2)
      iiniti(1) = 1
      iiniti(2) = 1
      ilasti(1) = k - 2
      ilasti(2) = k - 2
      dvm = mappl (plref, ahdr, axis, coeff, const, lvaddr, lvtype,
     +             iiniti, ilasti, istep, oiniti, olasti, ostep)
6     if (dopl (plref) .eq. 0) goto 7
      do 8 j = oiniti(1), olasti(1), ostep(1)
         do 8 i = oiniti(2), olasti(2), ostep(2)
            bptr(ahdr(3) + 1 + i + ahdr(2) * j) =
     +     ( bptr(ahdr(3) + 1 + (i - FDL2) + ahdr(2) * j) +
     +       bptr(ahdr(3) + 1 + i + ahdr(2) * (j - FDL1)) +
     +       bptr(ahdr(3) + 1 + (i + ADL2) + ahdr(2) * j) +
     +       bptr(ahdr(3) + 1 + i + ahdr(2) * (j + ADL1)) ) / 4.
8     continue
      goto 6
7     dvm = endpl (plref)
C     Sending exported elements
C     to support data dependence
      dvm = sendsh (shgref)
      dvm = waitsh (shgref)
C     Run-Time System termination
      dvm = lexit (0)
      end

PROGRAM IN C LANGUAGE

#include "dvmlib.h"
int main(int argc, char *argv[])
{ 
   long            InitPar = 0, UserRes = 0;
   long            StaticMV = 0, StaticDA = 0, ReDistrDA = 0,
                   StaticShG = 0, FullShd = 0, MaxShdCount = 1;
   long            ShdSignArray[2];
   long            i, j, Rank2 = 2, TypeSize = sizeof(float);
   AMRef           amref;
   PSRef           psref;
   AMViewRef       mvref;
   LoopRef         plref;
   ShadowGroupRef  fshgref, ashgref;
   AddrType        lvaddr[2];
   long            lvtype[2] = {0,0};
   long            amdim[2], disaxs[2], dispar[2];
   long            lshwdth[2], hshwdth[2], axis[2], coeff[2], cnst[2];
   long            iiniti[2], ilasti[2], istep[2];
   long            oiniti[2], olasti[2], ostep[2];
   long            ExtHdrSign = 0;
   DVMFILE         *OutFile;
   long k = 200;      /* size of each array dimension
                         with regular data dependence */
   long FDL1 = 4;     /* length of data
                         dependence for dimension 1 */
   long  FDL2 = 3;    /* length of data
                         dependence for dimension 2 */
   long ADL1 = 2;     /* length of data anti-dependence
                         for dimension 1 */
   long ADL2 = 1;     /* length of data anti-dependence
                         for dimension 2 */
   int f_inssh = 0;   /* flag to use inssh_ function to registrate
                         array in the shadow edge group if f_inssh = 0,
                         insshd_ function is used */
   long Ahdr[3]; /* header of array with
                    regular data dependence */
   rtl_init(InitPar, argc, argv); /* Run-Time System initialization */
   /* Creation of abstract machine representation
      and its mapping on processor subsystem */
   amref = getam_ ();
   psref = getps_ (&amref);
   amdim[0] = k;
   amdim[1] = k;
   mvref = crtamv_ (&amref, &Rank2, amdim, &StaticMV);
   disaxs[0] = 1;
   disaxs[1] = 2;
   dispar[0] = 0;
   dispar[1] = 0;
   distr_ (&mvref, &psref, &Rank2, disaxs, dispar);
   /* Creation and distribution of array
      with regular data dependence */
   lshwdth[0] = FDL1;   /* width of low shadow edge for dimension 1 is
                           length of data dependence for dimension 1 */
   lshwdth[1] = FDL2; /*   width of low shadow edge for dimension 2 is
                           length of data dependence for dimension 2 */
   hshwdth[0] = ADL1;    /* width of high shadow edge for
                            dimension 1 is length of data
                            anti-dependence for dimension 1 */
   hshwdth[1] = ADL2;    /* width of high shadow edge for
                            dimension 2 is length of data
                            anti-dependence for dimension 2 */
   crtda_ (Ahdr, &ExtHdrSign, NULL, &Rank2, &TypeSize, amdim,
           &StaticDA, &ReDistrDA, lshwdth, hshwdth);
   axis[0] = 1;
   axis[1] = 2;
   coeff[0] = 1;
   coeff[1] = 1;
   cnst[0] = 0;
   cnst[1] = 0;
   align_ (Ahdr,&mvref,axis,coeff,cnst);
   /*  Parallel loop of array initialization
       with regular data dependence */
   plref = crtpl_ (&Rank2);
   lvaddr[0] = (AddrType)&j;
   lvaddr[1] = (AddrType)&i;
   iiniti[0] = 0;
   iiniti[1] = 0;
   ilasti[0] = k - 1;
   ilasti[1] = k - 1;
   istep[0] = 1;
   istep[1] = 1;
   mappl_ (&plref, (PatternRef *)Ahdr, axis, coeff, cnst, lvaddr,
           lvtype, iiniti, ilasti, istep, oiniti, olasti, ostep);
   while (dopl_ (&plref))
      for ( j = oiniti[0]; j <= olasti[0]; j += ostep[0] )
         for ( i = oiniti[1]; i <= olasti[1]; i += ostep[1] )
             DAElm2(Ahdr, float, j, i) = 1.f;
   endpl_ (&plref);
   /*  Creation of shadow edge group
       to support data anti-dependence */
   ashgref = crtshg_ (&StaticShG);
   lshwdth[0] = 0;    /* zero width of low shadow edge for dimension 1 */
   lshwdth[1] = 0;    /* zero width of low shadow edge for dimension 2 */
   hshwdth[0] = ADL1; /* width of high shadow edge for
                         dimension 1 is length of data anti-dependence
                         for dimension 1 */
   hshwdth[1] = ADL2; /* width of high shadow edge
                         for dimension 2 is length
                         of data anti-dependence for dimension 2 */
   inssh_ (&ashgref, Ahdr, lshwdth, hshwdth, &FullShd);
    /* Initialization of shadow edges to support
      data anti-dependence */
   strtsh_ (&ashgref);
   waitsh_ (&ashgref);
   /*  Creation of shadow edge group
       to support data dependence */
   fshgref = crtshg_ (&StaticShG);
   lshwdth[0] = FDL1; /* width of low shadow edge for
                         dimension 1 is length of data
                         dependence for dimension 1 */
   lshwdth[1] = FDL2; /* width of low shadow edge for
                         dimension 2 is length of data
                         dependence for dimension 2 */
   hshwdth[0] = 0;    /* zero width of high shadow
                         edge for dimension 1 */
   hshwdth[1] = 0;    /* zero width of high shadow
                         edge for dimension 2 */
   if(f_inssh) /* inssh_ or insshd_ function is used */
   {
    /* inssh_ is used to registrate array in shadow edge group */
      inssh_ (&fshgref, Ahdr, lshwdth, hshwdth, &FullShd);
   }
   else
   {
     /* insshd_ is used to registrate array in shadow edge group */
     ShdSignArray[0] = 3; /* only low shadow edge is
                             used for dimension 1 */
     ShdSignArray[1] = 3; /* only low shadow edge is
                           used for dimension 2 */
     insshd_ (&fshgref, Ahdr, lshwdth, hshwdth, &MaxShdCount,
            ShdSignArray);
   }
   /* Receiving imported elements to support data dependence */
   recvsh_ (&fshgref);
   waitsh_ (&fshgref);
   /* PARALLEL LOOP WITH DATA DEPENDENCE AND
      ANTI-DEPENDENCE BETWEEN ITERATIONS */
   plref = crtpl_ (&Rank2);
   iiniti[0] = 1;
   iiniti[1] = 1;
   ilasti[0] = k - 2;
   ilasti[1] = k - 2;
   mappl_ (&plref, (PatternRef *)Ahdr, axis, coeff, cnst, lvaddr,
           lvtype, iiniti, ilasti, istep, oiniti, olasti, ostep);
   while (dopl_ (&plref))
      for ( j = oiniti[0]; j <= olasti[0]; j += ostep[0] )
         for ( i = oiniti[1]; i <= olasti[1]; i += ostep[1] )
               DAElm2(Ahdr, float, j, i) =
               ( DAElm2(Ahdr, float, j, i-FDL2) +
                 DAElm2(Ahdr, float, j-FDL1, i) +
                 DAElm2(Ahdr, float, j, i+ADL2) +
                 DAElm2(Ahdr, float, j+ADL1, i) ) / 4.f;
   endpl_ (&plref);
   /* Sending exported elements to support data dependence */
   sendsh_ (&fshgref);
   waitsh_ (&fshgref);
   /* Writing result array into the file
      and Run-Time System termination */
   OutFile = dvm_fopen("across.dat", "w+b");
   if(OutFile)
   { dvm_dfwrite(Ahdr, 0, OutFile);
     dvm_fclose(OutFile);
   }
   lexit_ (&UserRes);
   return (int)UserRes;
}

Lib-DVM interface description (contents) Part 1
(1-5)
Part 2
(6-7)
Part 3
(8-11)
Part 4
(12-13)
Part 5
(14-15)
Part 6
(16-18)
Part 7
(19)