#include <GlobalArray.h>
Public Methods | |
| GlobalArray (int type, int ndim, int dims[], char *arrayname, int chunk[]) | |
| Creates an ndim-dimensional array using the regular distribution model and returns integer handle representing the array. More... | |
| GlobalArray (int type, int ndim, int dims[], char *arrayname, int block[], int maps[]) | |
| Creates an array by following the user-specified distribution and returns integer handle representing the array. More... | |
| GlobalArray (int type, int ndim, int dims[], int width[], char *arrayname, int chunk[], char ghosts) | |
| Creates an ndim-dimensional array with a layer of ghost cells around the visible data on each processor using the regular distribution model and returns an integer handle representing the array. More... | |
| GlobalArray (int type, int ndim, int dims[], int width[], char *arrayname, int block[], int maps[], char ghosts) | |
| Creates an array with ghost cells by following the user-specified distribution and returns integer handle representing the array. More... | |
| GlobalArray (const GlobalArray &g_a, char *arrayname) | |
| Creates a new array by applying all the properties of another existing array. More... | |
| GlobalArray (const GlobalArray &g_a) | |
| Creates a new array by applying all the properties of another existing array. More... | |
| GlobalArray () | |
| Creates a 10x10 array of type "double"(default). More... | |
| ~GlobalArray () | |
| Destructor. More... | |
| int | handle () const |
| void | acc (int lo[], int hi[], void *buf, int ld[], void *alpha) const |
| Combines data from local array buffer with data in the global array section. More... | |
| void | access (int lo[], int hi[], void *ptr, int ld[]) const |
| Provides access to the specified patch of a global array. More... | |
| void | accessGhosts (int dims[], void *ptr, int ld[]) const |
| Provides access to the local patch of the global array. More... | |
| void | accessGhostElement (void *ptr, int subscript[], int ld[]) const |
| void | add (void *alpha, const GlobalArray *g_a, void *beta, const GlobalArray *g_b) const |
| The arrays are aded together elemet-wise: [for example: g_c.add(...,g_a, .., g_b);] c = alpha * a + beta * b The result c may replace one of he input arrays(a/b). More... | |
| void | addPatch (void *alpha, const GlobalArray *g_a, int alo[], int ahi[], void *beta, const GlobalArray *g_b, int blo[], int bhi[], int clo[], int chi[]) const |
| Patches of arrays (which must have the same number of elements) are added together element-wise. More... | |
| void | checkHandle (char *string) const |
| int | compareDistr (const GlobalArray *g_a) const |
| Compares distributions of two global arrays. More... | |
| void | copy (const GlobalArray *g_a) const |
| Copies elements in array represented by g_a into the array represented by g_b [say for example: g_b.copy(g_a);]. More... | |
| void | copyPatch (char trans, const GlobalArray *ga, int alo[], int ahi[], int blo[], int bhi[]) const |
| Copies elements in a patch of one array (ga) into another one (say for example:gb.copyPatch(...,ga,....); ). More... | |
| double | ddot (const GlobalArray *g_a) const |
| Computes element-wise dot product of the two arrays which must be of the same types and same number of elements. More... | |
| double | ddotPatch (char ta, int alo[], int ahi[], const GlobalArray *g_a, char tb, int blo[], int bhi[]) const |
| Computes the element-wise dot product of the two (possibly transposed) patches which must be of the same type and have the same number of elements. More... | |
| void | destroy () const |
| Deallocates the array and frees any associated resources. More... | |
| void | dgemm (char ta, char tb, int m, int n, int k, double alpha, const GlobalArray *g_a, const GlobalArray *g_b, double beta) const |
| Performs one of the matrix-matrix operations: [say: g_c.dgemm(..., g_a, g_b,..);]. More... | |
| void | diag (const GlobalArray *g_s, GlobalArray *g_v, void *eval) const |
| void | diagReuse (int control, const GlobalArray *g_s, GlobalArray *g_v, void *eval) const |
| Solve the generalized eigen-value problem returning all eigen-vectors and values in ascending order. More... | |
| void | diagStd (GlobalArray *g_v, void *eval) const |
| Solve the standard (non-generalized) eigenvalue problem returning all eigenvectors and values in the ascending order. More... | |
| void | diagSeq (const GlobalArray *g_s, const GlobalArray *g_v, void *eval) const |
| void | diagStdSeq (const GlobalArray *g_v, void *eval) const |
| void | distribution (int me, int *lo, int *hi) const |
| If no array elements are owned by process 'me', the range is returned as lo[]=0 and hi[]=-1 for all dimensions. More... | |
| float | fdot (const GlobalArray *g_a) const |
| float | fdotPatch (char t_a, int alo[], int ahi[], const GlobalArray *g_b, char t_b, int blo[], int bhi[]) const |
| void | fill (void *value) const |
| void | fillPatch (int lo[], int hi[], void *val) const |
| void | gather (void *v, int *subsarray[], int n) const |
| Gathers array elements from a global array into a local array. More... | |
| void | get (int lo[], int hi[], void *buf, int ld[]) const |
| One-side operations. More... | |
| int | hasGhosts () const |
| This function returns 1 if the global array has some dimensions for which the ghost cell width is greater than zero, it returns 0 otherwise. More... | |
| Integer | idot (const GlobalArray *g_a) const |
| Computes element-wise dot product of the two arrays which must be of the same types and same number of elements. More... | |
| long | idotPatch (char ta, int alo[], int ahi[], const GlobalArray *g_a, char tb, int blo[], int bhi[]) const |
| void | inquire (int *type, int *ndim, int dims[]) const |
| Returns data type and dimensions of the array. More... | |
| char * | inquireName () const |
| Returns the name of an array represented by the handle g_a. More... | |
| long | ldot (const GlobalArray *g_a) const |
| Computes element-wise dot product of the two arrays which must be of the same types and same number of elements. More... | |
| int | lltSolve (const GlobalArray *g_a) const |
| int | locate (int subscript[]) const |
| Return in owner the GA compute process id that 'owns' the data. More... | |
| int | locateRegion (int lo[], int hi[], int map[], int procs[]) const |
| Return the list of the GA processes id that 'own' the data. More... | |
| void | luSolve (char trans, const GlobalArray *g_a) const |
| void | matmulPatch (char transa, char transb, void *alpha, void *beta, const GlobalArray *g_a, int ailo, int aihi, int ajlo, int ajhi, const GlobalArray *g_b, int bilo, int bihi, int bjlo, int bjhi, int cilo, int cihi, int cjlo, int cjhi) const |
| void | matmulPatch (char transa, char transb, void *alpha, void *beta, const GlobalArray *g_a, int *alo, int *ahi, const GlobalArray *g_b, int *blo, int *bhi, int *clo, int *chi) const |
| N-dimensional Arrays:. More... | |
| void | nblock (int numblock[]) const |
| int | ndim () const |
| Returns the number of dimensions in array represented by the handle g_a. More... | |
| void | periodicAcc (int lo[], int hi[], void *buf, int ld[], void *alpha) const |
| void | periodicGet (int lo[], int hi[], void *buf, int ld[]) const |
| void | periodicPut (int lo[], int hi[], void *buf, int ld[]) const |
| void | print () const |
| Prints an entire array to the standard output. More... | |
| void | printDistribution () const |
| Prints the array distribution. More... | |
| void | printFile (FILE *file) const |
| Prints the array distribution to a file. More... | |
| void | printPatch (int *lo, int *hi, int pretty) const |
| Prints a patch of g_a array to the standard output. More... | |
| void | procTopology (int proc, int coord[]) const |
| void | put (int lo[], int hi[], void *buf, int ld[]) const |
| Copies data from local array buffer to the global array section . More... | |
| long | readInc (int subscript[], long inc) const |
| void | release (int lo[], int hi[]) const |
| void | releaseUpdate (int lo[], int hi[]) const |
| void | scale (void *value) const |
| Scales an array by the constant s. More... | |
| void | scalePatch (int lo[], int hi[], void *val) const |
| void | scatter (void *v, int *subsarray[], int n) const |
| Scatters array elements into a global array. More... | |
| void | sgemm (char ta, char tb, int m, int n, int k, float alpha, const GlobalArray *g_a, const GlobalArray *g_b, float beta) const |
| int | solve (const GlobalArray *g_a) const |
| int | spdInvert () const |
| It computes the inverse of a double precision using the Cholesky factorization of a NxN double precision symmetric positive definite matrix A stored in the global array represented by g_a. More... | |
| void | selectElem (char *op, void *val, int index[]) const |
| void | symmetrize () const |
| Symmmetrizes matrix A with handle A:=.5 * (A+A'). More... | |
| void | transpose (const GlobalArray *g_a) const |
| Transposes a matrix: B = A', where A and B are represented by handles g_a and g_b [say, g_b.transpose(g_a);]. More... | |
| void | updateGhosts () const |
| This call updates the ghost cell regions on each processor with the corresponding neighbor data from other processors. More... | |
| int | updateGhostDir (int dimension, int idir, int cflag) const |
| This function can be used to update the ghost cells along individual directions. More... | |
| DoubleComplex | zdot (const GlobalArray *g_a) const |
| Computes element-wise dot product of the two arrays which must be of the same types and same number of elements. More... | |
| DoubleComplex | zdotPatch (char ta, int alo[], int ahi[], const GlobalArray *g_a, char tb, int blo[], int bhi[]) const |
| void | zero () const |
| Sets value of all elements in the array to zero. More... | |
| void | zeroPatch (int lo[], int hi[]) const |
| void | zgemm (char ta, char tb, int m, int n, int k, DoubleComplex alpha, const GlobalArray *g_a, const GlobalArray *g_b, DoubleComplex beta) const |
| void | absValue () const |
| Take element-wise absolute value of the array. More... | |
| void | absValuePatch (int *lo, int *hi) const |
| Take element-wise absolute value of the patch. More... | |
| void | addConstant (void *alpha) const |
| Add the constant pointed by alpha to each element of the array. More... | |
| void | addConstantPatch (int *lo, int *hi, void *alpha) const |
| Add the constant pointed by alpha to each element of the patch. More... | |
| void | recip () const |
| Take element-wise reciprocal of the array. More... | |
| void | recipPatch (int *lo, int *hi) const |
| Take element-wise reciprocal of the patch. More... | |
| void | elemMultiply (const GlobalArray *g_a, const GlobalArray *g_b) const |
| Computes the element-wise product of the two arrays which must be of the same types and same number of elements. More... | |
| void | elemMultiplyPatch (const GlobalArray *g_a, int *alo, int *ahi, const GlobalArray *g_b, int *blo, int *bhi, int *clo, int *chi) const |
| Computes the element-wise product of the two patches which must be of the same types and same number of elements. More... | |
| void | elemDivide (const GlobalArray *g_a, const GlobalArray *g_b) const |
| void | elemDividePatch (const GlobalArray *g_a, int *alo, int *ahi, const GlobalArray *g_b, int *blo, int *bhi, int *clo, int *chi) const |
| Computes the element-wise quotient of the two patches which must be of the same types and same number of elements. More... | |
| void | elemMaximum (const GlobalArray *g_a, const GlobalArray *g_b) const |
| Computes the element-wise maximum of the two arrays which must be of the same types and same number of elements. More... | |
| void | elemMaximumPatch (const GlobalArray *g_a, int *alo, int *ahi, const GlobalArray *g_b, int *blo, int *bhi, int *clo, int *chi) const |
| Computes the element-wise maximum of the two patches which must be of the same types and same number of elements. More... | |
| void | elemMinimum (const GlobalArray *g_a, const GlobalArray *g_b) const |
| Computes the element-wise minimum of the two arrays which must be of the same types and same number of elements. More... | |
| void | elemMinimumPatch (const GlobalArray *g_a, int *alo, int *ahi, const GlobalArray *g_b, int *blo, int *bhi, int *clo, int *chi) const |
| Computes the element-wise minimum of the two patches which must be of the same types and same number of elements. More... | |
| void | stepMax (const GlobalArray *g_a, double *step) const |
| Calculates the largest multiple of a vector g_b that can be added to this vector g_a while keeping each element of this vector nonnegative. More... | |
| void | stepMax2 (const GlobalArray *g_vv, const GlobalArray *g_xxll, const GlobalArray *g_xxuu, double *step2) const |
| Calculates the largest step size that should be used in a projected bound line search. More... | |
| void | stepMaxPatch (int *alo, int *ahi, const GlobalArray *g_b, int *blo, int *bhi, double *step) const |
| void | stepMax2Patch (int *xxlo, int *xxhi, const GlobalArray *g_vv, int *vvlo, int *vvhi, const GlobalArray *g_xxll, int *xxlllo, int *xxllhi, const GlobalArray *g_xxuu, int *xxuulo, int *xxuuhi, double *step2) const |
| void | shiftDiagonal (void *c) const |
| Adds this constant to the diagonal elements of the matrix. More... | |
| void | setDiagonal (const GlobalArray *g_v) const |
| Sets the diagonal elements of this matrix g_a with the elements of the vector g_v. More... | |
| void | zeroDiagonal () const |
| Sets the diagonal elements of this matrix g_a with zeros. More... | |
| void | addDiagonal (const GlobalArray *g_v) const |
| Adds the elements of the vector g_v to the diagonal of this matrix g_a. More... | |
| void | getDiagonal (const GlobalArray *g_a) const |
| Inserts the diagonal elements of this matrix g_a into the vector g_v. More... | |
| void | scaleRows (const GlobalArray *g_v) const |
| Scales the rows of this matrix g_a using the vector g_v. More... | |
| void | scaleCols (const GlobalArray *g_v) const |
| Scales the columns of this matrix g_a using the vector g_v. More... | |
| void | norm1 (double *nm) const |
| Computes the 1-norm of the matrix or vector g_a. More... | |
| void | normInfinity (double *nm) const |
| Computes the 1-norm of the matrix or vector g_a. More... | |
| void | median (const GlobalArray *g_a, const GlobalArray *g_b, const GlobalArray *g_c) const |
| Computes the componentwise Median of three arrays g_a, g_b, and g_c, and stores the result in this array g_m. More... | |
| void | medianPatch (const GlobalArray *g_a, int *alo, int *ahi, const GlobalArray *g_b, int *blo, int *bhi, const GlobalArray *g_c, int *clo, int *chi, int *mlo, int *mhi) const |
| Computes the componentwise Median of three patches g_a, g_b, and g_c, and stores the result in this patch g_m. More... | |
| GlobalArray & | operator= (const GlobalArray &g_a) |
| int | operator== (const GlobalArray &g_a) const |
| int | operator!= (const GlobalArray &g_a) const |
Private Attributes | |
| int | mHandle |
Definition at line 7 of file GlobalArray.h.
|
||||||||||||||||||||||||
|
Creates an ndim-dimensional array using the regular distribution model and returns integer handle representing the array. The array can be distributed evenly or not. The control over the distribution is accomplished by specifying chunk (block) size for all or some of array dimensions. For example, for a 2-dimensional array, setting chunk[0]=dim[0] gives distribution by vertical strips (chunk[0]*dims[0]); setting chunk[1]=dim[1] gives distribution by horizontal strips (chunk[1]*dims[1]). Actual chunks will be modified so that they are at least the size of the minimum and each process has either zero or one chunk. Specifying chunk[i] as <1 will cause that dimension to be distributed evenly. As a convenience, when chunk is specified as NULL, the entire array is distributed evenly. This is a collective operation.
|
|
||||||||||||||||||||||||||||
|
Creates an array by following the user-specified distribution and returns integer handle representing the array. The distribution is specified as a Cartesian product of distributions for each dimension. The array indices start at 0. For example, the following figure demonstrates distribution of a 2-dimensional array 8x10 on 6 (or more) processors. nblock[2]={3,2}, the size of map array is s=5 and array map contains the following elements map={0,2,8, 0, 5}. The distribution is nonuniform because, P1 and P4 get 20 elements each and processors P0,P2,P3, and P5 only 10 elements each.
This is a collective operation.
|
|
||||||||||||||||||||||||||||||||
|
Creates an ndim-dimensional array with a layer of ghost cells around the visible data on each processor using the regular distribution model and returns an integer handle representing the array. The array can be distributed evenly or not evenly. The control over the distribution is accomplished by specifying chunk (block) size for all or some of the array dimensions. For example, for a 2-dimensional array, setting chunk(1)=dim(1) gives distribution by vertical strips (chunk(1)*dims(1)); setting chunk(2)=dim(2) gives distribution by horizontal strips (chunk(2)*dims(2)). Actual chunks will be modified so that they are at least the size of the minimum and each process has either zero or one chunk. Specifying chunk(i) as <1 will cause that dimension (i-th) to be distributed evenly. The width of the ghost cell layer in each dimension is specified using the array width(). The local data of the global array residing on each processor will have a layer width[n] ghosts cells wide on either side of the visible data along the dimension n.
|
|
||||||||||||||||||||||||||||||||||||
|
Creates an array with ghost cells by following the user-specified distribution and returns integer handle representing the array. The distribution is specified as a Cartesian product of distributions for each dimension. For example, the following figure demonstrates distribution of a 2-dimensional array 8x10 on 6 (or more) processors. nblock(2)={3,2}, the size of map array is s=5 and array map contains the following elements map={1,3,7, 1, 6}. The distribution is nonuniform because, P1 and P4 get 20 elements each and processors P0,P2,P3, and P5 only 10 elements each.
The array width[] is used to control the width of the ghost cell boundary around the visible data on each processor. The local data of the global array residing on each processor will have a layer width[n] ghosts cells wide on either side of the visible data along the dimension n. This is a collective operation.
|
|
||||||||||||
|
Creates a new array by applying all the properties of another existing array. This is a collective operation.
|
|
|
Creates a new array by applying all the properties of another existing array. This is a collective operation.
|
|
|
Creates a 10x10 array of type "double"(default).
|
|
|
Destructor.
|
|
|
Definition at line 172 of file GlobalArray.h. References mHandle.
00172 { return mHandle; }
|
|
||||||||||||||||||||||||
|
Combines data from local array buffer with data in the global array section. The local array is assumed to be have the same number of dimensions as the global array. global array section (lo[],hi[]) += *alpha * buffer This is a one-sided and atomic operation.
|
|
||||||||||||||||||||
|
Provides access to the specified patch of a global array. Returns array of leading dimensions ld and a pointer to the first element in the patch. This routine allows to access directly, in place elements in the local section of a global array. It useful for writing new GA operations. A call to ga_access normally follows a previous call to ga_distribution that returns coordinates of the patch associated with a processor. You need to make sure that the coordinates of the patch are valid (test values returned from ga_distribution). Each call to ga_access has to be followed by a call to either ga_release or ga_release_update. You can access in this fashion only local data. Since the data is shared with other processes, you need to consider issues of mutual exclusion. This operation is local.
|
|
||||||||||||||||
|
Provides access to the local patch of the global array. Returns leading dimension ld and and pointer for the data. This routine will provide access to the ghost cell data residing on each processor. Calls to NGA_Access_ghosts should normally follow a call to NGA_Distribution that returns coordinates of the visible data patch associated with a processor. You need to make sure that the coordinates of the patch are valid (test values returned from NGA_Distribution). You can only access local data. This is a local operation.
|
|
||||||||||||||||
|
|
|
||||||||||||||||||||
|
The arrays are aded together elemet-wise: [for example: g_c.add(...,g_a, .., g_b);] c = alpha * a + beta * b The result c may replace one of he input arrays(a/b). This is a collective operation. |
|
||||||||||||||||||||||||||||||||||||||||||||
|
Patches of arrays (which must have the same number of elements) are added together element-wise. c[ ][ ] = alpha * a[ ][ ] + beta * b[ ][ ]. This is a collective operation.
|
|
|
|
|
|
Compares distributions of two global arrays. Returns 0 if distributions are identical and 1 when they are not. This is a collective operation.
|
|
|
Copies elements in array represented by g_a into the array represented by g_b [say for example: g_b.copy(g_a);]. The arrays must be the same type, shape, and identically aligned. This is a collective operation.
|
|
||||||||||||||||||||||||||||
|
Copies elements in a patch of one array (ga) into another one (say for example:gb.copyPatch(...,ga,....); ). The patches of arrays may be of different shapes but must have the same number of elements. Patches must be nonoverlapping (if gb=ga). trans = 'N' or 'n' means that the transpose operator should not be applied. trans = 'T' or 't' means that transpose operator should be applied. This is a collective operation.
|
|
|
Computes element-wise dot product of the two arrays which must be of the same types and same number of elements. return value = SUM_ij a(i,j)*b(i,j) This is a collective operation.
|
|
||||||||||||||||||||||||||||||||
|
Computes the element-wise dot product of the two (possibly transposed) patches which must be of the same type and have the same number of elements.
|
|
|
Deallocates the array and frees any associated resources.
|
|
||||||||||||||||||||||||||||||||||||||||
|
Performs one of the matrix-matrix operations: [say: g_c.dgemm(..., g_a, g_b,..);].
C := alpha*op( A )*op( B ) + beta*C,
|
|
||||||||||||||||
|
|
|
||||||||||||||||||||
|
Solve the generalized eigen-value problem returning all eigen-vectors and values in ascending order. Recommended for REPEATED calls if g_s is unchanged. Values of the control flag: value action/purpose 0 indicates first call to the eigensolver >0 consecutive calls (reuses factored g_s) <0 only erases factorized g_s; g_v and eval unchanged (should be called after previous use if another eigenproblem, i.e., different g_a and g_s, is to be solved) The input matrices are not destroyed. This is a collective operation.
|
|
||||||||||||
|
Solve the standard (non-generalized) eigenvalue problem returning all eigenvectors and values in the ascending order. The input matrix is neither overwritten nor destroyed. This is a collective operation.
|
|
||||||||||||||||
|
|
|
||||||||||||
|
|
|
||||||||||||||||
|
If no array elements are owned by process 'me', the range is returned as lo[]=0 and hi[]=-1 for all dimensions. The operation is local.
|
|
|
|
|
||||||||||||||||||||||||||||||||
|
|
|
|
|
|
||||||||||||||||
|
|
|
||||||||||||||||
|
Gathers array elements from a global array into a local array. The contents of the input arrays (v, subscrArray) are preserved, but their contents might be (consistently) shuffled on return. for(k=0; k<= n; k++){ v[k] = a[subsArray[k][0]][subsArray[k][1]][subsArray[k][2]]...; } This is a one-sided operation.
|
|
||||||||||||||||||||
|
One-side operations. Copies data from global array section to the local array buffer. The local array is assumed to be have the same number of dimensions as the global array. Any detected inconsitencies/errors in the input arguments are fatal. Example: For ga_get operation transfering data from the [10:14,0:4] section of 2-dimensional 15x10 global array into local buffer 5x10 array we have: lo={10,0}, hi={14,4}, ld={10}
|
|
|
This function returns 1 if the global array has some dimensions for which the ghost cell width is greater than zero, it returns 0 otherwise. This is a collective operation. |
|
|
Computes element-wise dot product of the two arrays which must be of the same types and same number of elements. return value = SUM_ij a(i,j)*b(i,j) This is a collective operation.
|
|
||||||||||||||||||||||||||||||||
|
|
|
||||||||||||||||
|
Returns data type and dimensions of the array. This operation is local.
|
|
|
Returns the name of an array represented by the handle g_a. This operation is local. |
|
|
Computes element-wise dot product of the two arrays which must be of the same types and same number of elements. return value = SUM_ij a(i,j)*b(i,j) This is a collective operation.
|
|
|
= 0 : successful exit > 0 : the leading minor of this order is not positive definite and the factorization could not be completed This is a collective operation. |
|
|
Return in owner the GA compute process id that 'owns' the data. If any element of subscript[] is out of bounds "-1" is returned. This operation is local.
|
|
||||||||||||||||||||
|
Return the list of the GA processes id that 'own' the data. Parts of the specified patch might be actually 'owned' by several processes. If lo/hi are out of bounds "0" is returned, otherwise return value is equal to the number of processes that hold the data. This operation is local. map[i][0:ndim-1] - lo[i] map[i][ndim:2*ndim-1] - hi[i] procs[i] - processor id that owns data in patch lo[i]:hi[i]
|
|
||||||||||||
|
op(A) = A or A' depending on the parameter trans: trans = 'N' or 'n' means that the transpose operator should not be applied. trans = 'T' or 't' means that the transpose operator should be applied. Matrix A is a general real matrix. Matrix B contains possibly multiple rhs vectors. The array associated with the handle g_b is overwritten by the solution matrix X. This is a collective operation. |
|
||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||
|
C[cilo:cihi,cjlo:cjhi] := alpha* AA[ailo:aihi,ajlo:ajhi] * BB[bilo:bihi,bjlo:bjhi] ) + beta*C[cilo:cihi,cjlo:cjhi], where AA = op(A), BB = op(B), and op( X ) is one of op( X ) = X or op( X ) = X', Valid values for transpose arguments: 'n', 'N', 't', 'T'. It works for both double and DoubleComplex data tape. This is a collective operation. |
|
||||||||||||||||||||||||||||||||||||||||||||||||||||
|
N-dimensional Arrays:.
C[clo[]:chi[]] := alpha* AA[alo[]:ahi[]] * BB[blo[]:bhi[]]) + beta*C[clo[]:chi[]], where AA = op(A), BB = op(B), and op( X ) is one of op( X ) = X or op( X ) = X', Valid values for transpose arguments: 'n', 'N', 't', 'T'. It works for both double and DoubleComplex data tape. This is a collective operation. |
|
|
|
|
|
Returns the number of dimensions in array represented by the handle g_a. This operation is local. |
|
||||||||||||||||||||||||
|
|
|
||||||||||||||||||||
|
|
|
||||||||||||||||||||
|
|
|
|
Prints an entire array to the standard output. This is a collective operation. |
|
|
Prints the array distribution. This is a collective operation. |
|
|
Prints the array distribution to a file. This is a collective operation. |
|
||||||||||||||||
|
Prints a patch of g_a array to the standard output. If pretty has the value 0 then output is printed in a dense fashion. If pretty has the value 1 then output is formatted and rows/columns labeled. This is a collective operation.
|
|
||||||||||||
|
|
|
||||||||||||||||||||
|
Copies data from local array buffer to the global array section . The local array is assumed to be have the same number of dimensions as the global array. Any detected inconsitencies/errors in input arguments are fatal. This is a one-sided operation.
|
|
||||||||||||
|
*BEGIN CRITICAL SECTION* old_value = a(subscript) a(subscript) += inc *END CRITICAL SECTION* return old_value This is a one-sided and atomic operation. |
|
||||||||||||
|
NGA_Distribution(g_a, myproc, lo,hi); NGA_Access(g_a, lo, hi, &ptr, ld); <operate on the data referenced by ptr> GA_Release(g_a, lo, hi); NOTE: see restrictions specified for ga_access. This operation is local. |
|
||||||||||||
|
|
|
|
Scales an array by the constant s. Note that the library is unable to detect errors when the pointed value is of different type than the array. This is a collective operation.
|
|
||||||||||||||||
|
|
|
||||||||||||||||
|
Scatters array elements into a global array. The contents of the input arrays (v,subscrArray) are preserved, but their contents might be (consistently) shuffled on return. for(k=0; k<= n; k++){ a[subsArray[k][0]][subsArray[k][1]][subsArray[k][2]]... = v[k]; } This is a one-sided operation.
|
|
||||||||||||||||||||||||||||||||||||||||
|
C := alpha*op( A )*op( B ) + beta*C, where op( X ) is one of op( X ) = X or op( X ) = X', alpha and beta are scalars, and A, B and C are matrices, with op( A ) an m by k matrix, op( B ) a k by n matrix and C an m by n matrix. On entry, transa specifies the form of op( A ) to be used in the matrix multiplication as follows: ta = 'N' or 'n', op( A ) = A. ta = 'T' or 't', op( A ) = A'. This is a collective operation. |
|
|
= 0 : Cholesky factoriztion was succesful > 0 : the leading minor of this order is not positive definite, Cholesky factorization could not be completed and LU factoriztion was used This is a collective operation. |
|
|
It computes the inverse of a double precision using the Cholesky factorization of a NxN double precision symmetric positive definite matrix A stored in the global array represented by g_a. On successful exit, A will contain the inverse. It returns = 0 : successful exit * > 0 : the leading minor of this order is not positive definite and the factorization could not be completed < 0 : it returns the index i of the (i,i) element of the factor L/U that is zero and the inverse could not be computed This is a collective operation. |
|
||||||||||||||||
|
|
|
|
Symmmetrizes matrix A with handle A:=.5 * (A+A'). This is a collective operation |
|
|
Transposes a matrix: B = A', where A and B are represented by handles g_a and g_b [say, g_b.transpose(g_a);]. This is a collective operation. |
|
|
This call updates the ghost cell regions on each processor with the corresponding neighbor data from other processors. The operation assumes that all data is wrapped around using periodic boundary data so that ghost cell data that goes beyound an array boundary is wrapped around to the other end of the array. The GA_Update_ghosts call contains two GA_Sync calls before and after the actual update operation. For some applications these calls may be unecessary, if so they can be removed using the GA_Mask_sync subroutine. This is a collective operation. |
|
||||||||||||||||
|
This function can be used to update the ghost cells along individual directions. It is designed for algorithms that can overlap updates with computation. The variable dimension indicates which coordinate direction is to be updated (e.g. dimension = 1 would correspond to the y axis in a two or three dimensional system), the variable idir can take the values +/-1 and indicates whether the side that is to be updated lies in the positive or negative direction, and cflag indicates whether or not the corners on the side being updated are to be included in the update. The following calls would be equivalent to a call to GA_Update_ghosts for a 2-dimensional system:
status = NGA_Update_ghost_dir(g_a,0,-1,1); The variable cflag is set equal to 1 (or non-zero) in the first two calls so that the corner ghost cells are update, it is set equal to 0 in the second two calls to avoid redundant updates of the corners. Note that updating the ghosts cells using several independent calls to the nga_update_ghost_dir functions is generally not as efficient as using GA_Update_ghosts unless the individual calls can be effectively overlapped with computation. This is a collective operation.
|
|
|
Computes element-wise dot product of the two arrays which must be of the same types and same number of elements. return value = SUM_ij a(i,j)*b(i,j) This is a collective operation.
|
|
||||||||||||||||||||||||||||||||
|
|
|
|
Sets value of all elements in the array to zero. This is a collective operation. |
|
||||||||||||
|
|
|
||||||||||||||||||||||||||||||||||||||||
|
ta = 'N' or 'n', op( A ) = A. ta = 'T' or 't', op( A ) = A'. * This is a collective operation. |
|
|
Take element-wise absolute value of the array. This is a collective operation. |
|
||||||||||||
|
Take element-wise absolute value of the patch. This is a collective operation.
|
|
|
Add the constant pointed by alpha to each element of the array. This is a collective operation.
|
|
||||||||||||||||
|
Add the constant pointed by alpha to each element of the patch. This is a collective operation.
|
|
|
Take element-wise reciprocal of the array. This is a collective operation. |
|
||||||||||||
|
Take element-wise reciprocal of the patch. This is a collective operation.
|
|
||||||||||||
|
Computes the element-wise product of the two arrays which must be of the same types and same number of elements. For two-dimensional arrays, c(i, j) = a(i,j)*b(i,j) The result (c) may replace one of the input arrays (a/b). This is a collective operation.
|
|
||||||||||||||||||||||||||||||||||||
|
Computes the element-wise product of the two patches which must be of the same types and same number of elements. For two-dimensional arrays, c(i, j) = a(i,j)*b(i,j) The result (c) may replace one of the input arrays (a/b). This is a collective operation.
|
|
||||||||||||
|
c(i, j) = a(i,j)/b(i,j) The result (c) may replace one of the input arrays (a/b). If one of the elements of array g_b is zero, the quotient for the element of g_c will be set to GA_NEGATIVE_INFINITY. This is a collective operation. |
|
||||||||||||||||||||||||||||||||||||
|
Computes the element-wise quotient of the two patches which must be of the same types and same number of elements. For two-dimensional arrays, c(i, j) = a(i,j)/b(i,j) The result (c) may replace one of the input arrays (a/b). This is a collective operation.
|
|
||||||||||||
|
Computes the element-wise maximum of the two arrays which must be of the same types and same number of elements. For two dimensional arrays, c(i, j) = max{a(i,j), b(i,j)} The result (c) may replace one of the input arrays (a/b). This is a collective operation.
|
|
||||||||||||||||||||||||||||||||||||
|
Computes the element-wise maximum of the two patches which must be of the same types and same number of elements. For two-dimensional of noncomplex arrays, c(i, j) = max{a(i,j), b(i,j)} If the data type is complex, then c(i, j).real = max{ |a(i,j)|, |b(i,j)|} while c(i,j).image = 0. The result (c) may replace one of the input arrays (a/b). This is a collective operation.
|
|
||||||||||||
|
Computes the element-wise minimum of the two arrays which must be of the same types and same number of elements. For two dimensional arrays, c(i, j) = min{a(i,j), b(i,j)} The result (c) may replace one of the input arrays (a/b). This is a collective operation.
|
|
||||||||||||||||||||||||||||||||||||
|
Computes the element-wise minimum of the two patches which must be of the same types and same number of elements. For two-dimensional of noncomplex arrays, c(i, j) = min{a(i,j), b(i,j)} If the data type is complex, then c(i, j).real = min{ |a(i,j)|, |b(i,j)|} while c(i,j).image = 0. The result (c) may replace one of the input arrays (a/b). This is a collective operation.
|
|
||||||||||||
|
Calculates the largest multiple of a vector g_b that can be added to this vector g_a while keeping each element of this vector nonnegative. This is a collective operation.
|
|
||||||||||||||||||||
|
Calculates the largest step size that should be used in a projected bound line search. This is a collective operation.
|
|
||||||||||||||||||||||||||||
|
|
|
||||||||||||||||||||||||||||||||||||||||||||||||||||
|
|
|
|
Adds this constant to the diagonal elements of the matrix. This is a collective operation.
|
|
|
Sets the diagonal elements of this matrix g_a with the elements of the vector g_v. This is a collective operation.
|
|
|
Sets the diagonal elements of this matrix g_a with zeros. This is a collective operation. |
|
|
Adds the elements of the vector g_v to the diagonal of this matrix g_a. This is a collective operation.
|
|
|
Inserts the diagonal elements of this matrix g_a into the vector g_v. This is a collective operation.
|
|
|
Scales the rows of this matrix g_a using the vector g_v. This is a collective operation.
|
|
|
Scales the columns of this matrix g_a using the vector g_v. This is a collective operation.
|
|
|
Computes the 1-norm of the matrix or vector g_a. This is a collective operation.
|
|
|
Computes the 1-norm of the matrix or vector g_a. This is a collective operation.
|
|
||||||||||||||||
|
Computes the componentwise Median of three arrays g_a, g_b, and g_c, and stores the result in this array g_m. The result (m) may replace one of the input arrays (a/b/c). This is a collective operation.
|
|
||||||||||||||||||||||||||||||||||||||||||||||||
|
Computes the componentwise Median of three patches g_a, g_b, and g_c, and stores the result in this patch g_m. The result (m) may replace one of the input patches (a/b/c). This is a collective operation.
|
|
|
|
|
|
|
|
|
|
|
|
Definition at line 1403 of file GlobalArray.h. Referenced by handle. |
1.2.14 written by Dimitri van Heesch,
© 1997-2002