Actual source code: mpicusp.cu

petsc-3.6.4 2016-04-12
Report Typos and Errors
  2: /*
  3:    This file contains routines for Parallel vector operations.
  4:  */
  5: #define PETSC_SKIP_COMPLEX

  7: #include <petscconf.h>
  8: PETSC_CUDA_EXTERN_C_BEGIN
  9: #include <../src/vec/vec/impls/mpi/pvecimpl.h>   /*I  "petscvec.h"   I*/
 10: PETSC_CUDA_EXTERN_C_END
 11: #include <../src/vec/vec/impls/seq/seqcusp/cuspvecimpl.h>

 15: PetscErrorCode VecDestroy_MPICUSP(Vec v)
 16: {

 20:   try {
 21:     if (v->spptr) {
 22:       delete ((Vec_CUSP*)v->spptr)->GPUarray;
 23:       delete (Vec_CUSP*) v->spptr;
 24:     }
 25:   } catch(char * ex) {
 26:     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CUSP error: %s", ex);
 27:   }
 28:   VecDestroy_MPI(v);
 29:   return(0);
 30: }

 34: PetscErrorCode VecNorm_MPICUSP(Vec xin,NormType type,PetscReal *z)
 35: {
 36:   PetscReal      sum,work = 0.0;

 40:   if (type == NORM_2 || type == NORM_FROBENIUS) {
 41:     VecNorm_SeqCUSP(xin,NORM_2,&work);
 42:     work *= work;
 43:     MPI_Allreduce(&work,&sum,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)xin));
 44:     *z    = PetscSqrtReal(sum);
 45:     //printf("VecNorm_MPICUSP : z=%1.5g\n",*z);
 46:   } else if (type == NORM_1) {
 47:     /* Find the local part */
 48:     VecNorm_SeqCUSP(xin,NORM_1,&work);
 49:     /* Find the global max */
 50:     MPI_Allreduce(&work,z,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)xin));
 51:   } else if (type == NORM_INFINITY) {
 52:     /* Find the local max */
 53:     VecNorm_SeqCUSP(xin,NORM_INFINITY,&work);
 54:     /* Find the global max */
 55:     MPI_Allreduce(&work,z,1,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)xin));
 56:   } else if (type == NORM_1_AND_2) {
 57:     PetscReal temp[2];
 58:     VecNorm_SeqCUSP(xin,NORM_1,temp);
 59:     VecNorm_SeqCUSP(xin,NORM_2,temp+1);
 60:     temp[1] = temp[1]*temp[1];
 61:     MPI_Allreduce(temp,z,2,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)xin));
 62:     z[1] = PetscSqrtReal(z[1]);
 63:   }
 64:   return(0);
 65: }

 69: PetscErrorCode VecDot_MPICUSP(Vec xin,Vec yin,PetscScalar *z)
 70: {
 71:   PetscScalar    sum,work;

 75:   VecDot_SeqCUSP(xin,yin,&work);
 76:   MPI_Allreduce(&work,&sum,1,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)xin));
 77:   *z   = sum;
 78:   return(0);
 79: }

 83: PetscErrorCode VecTDot_MPICUSP(Vec xin,Vec yin,PetscScalar *z)
 84: {
 85:   PetscScalar    sum,work;

 89:   VecTDot_SeqCUSP(xin,yin,&work);
 90:   MPI_Allreduce(&work,&sum,1,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)xin));
 91:   *z   = sum;
 92:   return(0);
 93: }

 97: PetscErrorCode VecMDot_MPICUSP(Vec xin,PetscInt nv,const Vec y[],PetscScalar *z)
 98: {
 99:   PetscScalar    awork[128],*work = awork;

103:   if (nv > 128) {
104:     PetscMalloc1(nv,&work);
105:   }
106:   VecMDot_SeqCUSP(xin,nv,y,work);
107:   MPI_Allreduce(work,z,nv,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)xin));
108:   if (nv > 128) {
109:     PetscFree(work);
110:   }
111:   return(0);
112: }

114: /*MC
115:    VECMPICUSP - VECMPICUSP = "mpicusp" - The basic parallel vector, modified to use CUSP

117:    Options Database Keys:
118: . -vec_type mpicusp - sets the vector type to VECMPICUSP during a call to VecSetFromOptions()

120:   Level: beginner

122: .seealso: VecCreate(), VecSetType(), VecSetFromOptions(), VecCreateMpiWithArray(), VECMPI, VecType, VecCreateMPI(), VecCreateMpi()
123: M*/


128: PetscErrorCode VecDuplicate_MPICUSP(Vec win,Vec *v)
129: {
131:   Vec_MPI        *vw,*w = (Vec_MPI*)win->data;
132:   PetscScalar    *array;

135:   VecCreate(PetscObjectComm((PetscObject)win),v);
136:   PetscLayoutReference(win->map,&(*v)->map);

138:   VecCreate_MPI_Private(*v,PETSC_FALSE,w->nghost,0);
139:   vw   = (Vec_MPI*)(*v)->data;
140:   PetscMemcpy((*v)->ops,win->ops,sizeof(struct _VecOps));

142:   /* save local representation of the parallel vector (and scatter) if it exists */
143:   if (w->localrep) {
144:     VecGetArray(*v,&array);
145:     VecCreateSeqWithArray(PETSC_COMM_SELF,1,win->map->n+w->nghost,array,&vw->localrep);
146:     PetscMemcpy(vw->localrep->ops,w->localrep->ops,sizeof(struct _VecOps));
147:     VecRestoreArray(*v,&array);
148:     PetscLogObjectParent((PetscObject)*v,(PetscObject)vw->localrep);
149:     vw->localupdate = w->localupdate;
150:     if (vw->localupdate) {
151:       PetscObjectReference((PetscObject)vw->localupdate);
152:     }
153:   }

155:   /* New vector should inherit stashing property of parent */
156:   (*v)->stash.donotstash   = win->stash.donotstash;
157:   (*v)->stash.ignorenegidx = win->stash.ignorenegidx;

159:   /* change type_name appropriately */
160:   PetscObjectChangeTypeName((PetscObject)(*v),VECMPICUSP);

162:   PetscObjectListDuplicate(((PetscObject)win)->olist,&((PetscObject)(*v))->olist);
163:   PetscFunctionListDuplicate(((PetscObject)win)->qlist,&((PetscObject)(*v))->qlist);
164:   (*v)->map->bs   = PetscAbs(win->map->bs);
165:   (*v)->bstash.bs = win->bstash.bs;
166:   return(0);
167: }

171: PetscErrorCode VecDotNorm2_MPICUSP(Vec s,Vec t,PetscScalar *dp,PetscScalar *nm)
172: {
174:   PetscScalar    work[2],sum[2];

177:   VecDotNorm2_SeqCUSP(s,t,work,work+1);
178:   MPI_Allreduce(&work,&sum,2,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)s));
179:   *dp  = sum[0];
180:   *nm  = sum[1];
181:   //printf("VecDotNorm2_MPICUSP=%1.5g,%1.5g\n",PetscRealPart(*dp),PetscImaginaryPart(*dp));
182:   //printf("VecDotNorm2_MPICUSP=%1.5g,%1.5g\n",PetscRealPart(*nm),PetscImaginaryPart(*nm));
183:   return(0);
184: }

188: PETSC_EXTERN PetscErrorCode VecCreate_MPICUSP(Vec vv)
189: {

193:   VecCreate_MPI_Private(vv,PETSC_FALSE,0,0);
194:   PetscObjectChangeTypeName((PetscObject)vv,VECMPICUSP);

196:   vv->ops->dotnorm2               = VecDotNorm2_MPICUSP;
197:   vv->ops->waxpy                  = VecWAXPY_SeqCUSP;
198:   vv->ops->duplicate              = VecDuplicate_MPICUSP;
199:   vv->ops->dot                    = VecDot_MPICUSP;
200:   vv->ops->mdot                   = VecMDot_MPICUSP;
201:   vv->ops->tdot                   = VecTDot_MPICUSP;
202:   vv->ops->norm                   = VecNorm_MPICUSP;
203:   vv->ops->scale                  = VecScale_SeqCUSP;
204:   vv->ops->copy                   = VecCopy_SeqCUSP;
205:   vv->ops->set                    = VecSet_SeqCUSP;
206:   vv->ops->swap                   = VecSwap_SeqCUSP;
207:   vv->ops->axpy                   = VecAXPY_SeqCUSP;
208:   vv->ops->axpby                  = VecAXPBY_SeqCUSP;
209:   vv->ops->maxpy                  = VecMAXPY_SeqCUSP;
210:   vv->ops->aypx                   = VecAYPX_SeqCUSP;
211:   vv->ops->axpbypcz               = VecAXPBYPCZ_SeqCUSP;
212:   vv->ops->pointwisemult          = VecPointwiseMult_SeqCUSP;
213:   vv->ops->setrandom              = VecSetRandom_SeqCUSP;
214:   vv->ops->replacearray           = VecReplaceArray_SeqCUSP;
215:   vv->ops->dot_local              = VecDot_SeqCUSP;
216:   vv->ops->tdot_local             = VecTDot_SeqCUSP;
217:   vv->ops->norm_local             = VecNorm_SeqCUSP;
218:   vv->ops->mdot_local             = VecMDot_SeqCUSP;
219:   vv->ops->destroy                = VecDestroy_MPICUSP;
220:   vv->ops->pointwisedivide        = VecPointwiseDivide_SeqCUSP;
221:   vv->ops->getlocalvector         = VecGetLocalVector_SeqCUSP;
222:   vv->ops->restorelocalvector     = VecRestoreLocalVector_SeqCUSP;
223:   vv->ops->getlocalvectorread     = VecGetLocalVector_SeqCUSP;
224:   vv->ops->restorelocalvectorread = VecRestoreLocalVector_SeqCUSP;
225:   /* place array?
226:      reset array?
227:      get values?
228:   */
229:   VecCUSPAllocateCheck(vv);CHKERRCUSP(ierr);
230:   vv->valid_GPU_array      = PETSC_CUSP_GPU;
231:   VecSet(vv,0.0);
232:   return(0);
233: }

237: PETSC_EXTERN PetscErrorCode VecCreate_CUSP(Vec v)
238: {
240:   PetscMPIInt    size;

243:   MPI_Comm_size(PetscObjectComm((PetscObject)v),&size);
244:   if (size == 1) {
245:     VecSetType(v,VECSEQCUSP);
246:   } else {
247:     VecSetType(v,VECMPICUSP);
248:   }
249:   return(0);
250: }