Actual source code: mpiviennacl.cxx

petsc-3.7.1 2016-05-15
Report Typos and Errors
  2: /*
  3:    This file contains routines for Parallel vector operations.
  4:  */
  5: #include <petscconf.h>
  6: #include <../src/vec/vec/impls/mpi/pvecimpl.h>   /*I  "petscvec.h"   I*/
  7: #include <../src/vec/vec/impls/seq/seqviennacl/viennaclvecimpl.h>

 11: PetscErrorCode VecDestroy_MPIViennaCL(Vec v)
 12: {

 16:   try {
 17:     if (v->spptr) {
 18:       delete ((Vec_ViennaCL*)v->spptr)->GPUarray;
 19:       delete (Vec_ViennaCL*) v->spptr;
 20:     }
 21:   } catch(std::exception const & ex) {
 22:     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"ViennaCL error: %s", ex.what());
 23:   }
 24:   VecDestroy_MPI(v);
 25:   return(0);
 26: }

 30: PetscErrorCode VecNorm_MPIViennaCL(Vec xin,NormType type,PetscReal *z)
 31: {
 32:   PetscReal      sum,work = 0.0;

 36:   if (type == NORM_2 || type == NORM_FROBENIUS) {
 37:     VecNorm_SeqViennaCL(xin,NORM_2,&work);
 38:     work *= work;
 39:     MPIU_Allreduce(&work,&sum,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)xin));
 40:     *z    = PetscSqrtReal(sum);
 41:   } else if (type == NORM_1) {
 42:     /* Find the local part */
 43:     VecNorm_SeqViennaCL(xin,NORM_1,&work);
 44:     /* Find the global max */
 45:     MPIU_Allreduce(&work,z,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)xin));
 46:   } else if (type == NORM_INFINITY) {
 47:     /* Find the local max */
 48:     VecNorm_SeqViennaCL(xin,NORM_INFINITY,&work);
 49:     /* Find the global max */
 50:     MPIU_Allreduce(&work,z,1,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)xin));
 51:   } else if (type == NORM_1_AND_2) {
 52:     PetscReal temp[2];
 53:     VecNorm_SeqViennaCL(xin,NORM_1,temp);
 54:     VecNorm_SeqViennaCL(xin,NORM_2,temp+1);
 55:     temp[1] = temp[1]*temp[1];
 56:     MPIU_Allreduce(temp,z,2,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)xin));
 57:     z[1] = PetscSqrtReal(z[1]);
 58:   }
 59:   return(0);
 60: }

 64: PetscErrorCode VecDot_MPIViennaCL(Vec xin,Vec yin,PetscScalar *z)
 65: {
 66:   PetscScalar    sum,work;

 70:   VecDot_SeqViennaCL(xin,yin,&work);
 71:   MPIU_Allreduce(&work,&sum,1,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)xin));
 72:   *z   = sum;
 73:   return(0);
 74: }

 78: PetscErrorCode VecTDot_MPIViennaCL(Vec xin,Vec yin,PetscScalar *z)
 79: {
 80:   PetscScalar    sum,work;

 84:   VecTDot_SeqViennaCL(xin,yin,&work);
 85:   MPIU_Allreduce(&work,&sum,1,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)xin));
 86:   *z   = sum;
 87:   return(0);
 88: }

 92: PetscErrorCode VecMDot_MPIViennaCL(Vec xin,PetscInt nv,const Vec y[],PetscScalar *z)
 93: {
 94:   PetscScalar    awork[128],*work = awork;

 98:   if (nv > 128) {
 99:     PetscMalloc1(nv,&work);
100:   }
101:   VecMDot_SeqViennaCL(xin,nv,y,work);
102:   MPIU_Allreduce(work,z,nv,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)xin));
103:   if (nv > 128) {
104:     PetscFree(work);
105:   }
106:   return(0);
107: }

109: /*MC
110:    VECMPIVIENNACL - VECMPIVIENNACL = "mpiviennacl" - The basic parallel vector, modified to use ViennaCL

112:    Options Database Keys:
113: . -vec_type mpiviennacl - sets the vector type to VECMPIVIENNACL during a call to VecSetFromOptions()

115:   Level: beginner

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


123: PetscErrorCode VecDuplicate_MPIViennaCL(Vec win,Vec *v)
124: {
126:   Vec_MPI        *vw,*w = (Vec_MPI*)win->data;
127:   PetscScalar    *array;

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

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

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

150:   /* New vector should inherit stashing property of parent */
151:   (*v)->stash.donotstash   = win->stash.donotstash;
152:   (*v)->stash.ignorenegidx = win->stash.ignorenegidx;

154:   /* change type_name appropriately */
155:   PetscObjectChangeTypeName((PetscObject)(*v),VECMPIVIENNACL);

157:   PetscObjectListDuplicate(((PetscObject)win)->olist,&((PetscObject)(*v))->olist);
158:   PetscFunctionListDuplicate(((PetscObject)win)->qlist,&((PetscObject)(*v))->qlist);
159:   (*v)->map->bs   = PetscAbs(win->map->bs);
160:   (*v)->bstash.bs = win->bstash.bs;
161:   return(0);
162: }

166: PetscErrorCode VecDotNorm2_MPIViennaCL(Vec s,Vec t,PetscScalar *dp,PetscScalar *nm)
167: {
169:   PetscScalar    work[2],sum[2];

172:   VecDotNorm2_SeqViennaCL(s,t,work,work+1);
173:   MPIU_Allreduce((void*)&work,(void*)&sum,2,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)s));
174:   *dp  = sum[0];
175:   *nm  = sum[1];
176:   return(0);
177: }

181: PETSC_EXTERN PetscErrorCode VecCreate_MPIViennaCL(Vec vv)
182: {

186:   VecCreate_MPI_Private(vv,PETSC_FALSE,0,0);
187:   PetscObjectChangeTypeName((PetscObject)vv,VECMPIVIENNACL);

189:   vv->ops->dotnorm2        = VecDotNorm2_MPIViennaCL;
190:   vv->ops->waxpy           = VecWAXPY_SeqViennaCL;
191:   vv->ops->duplicate       = VecDuplicate_MPIViennaCL;
192:   vv->ops->dot             = VecDot_MPIViennaCL;
193:   vv->ops->mdot            = VecMDot_MPIViennaCL;
194:   vv->ops->tdot            = VecTDot_MPIViennaCL;
195:   vv->ops->norm            = VecNorm_MPIViennaCL;
196:   vv->ops->scale           = VecScale_SeqViennaCL;
197:   vv->ops->copy            = VecCopy_SeqViennaCL;
198:   vv->ops->set             = VecSet_SeqViennaCL;
199:   vv->ops->swap            = VecSwap_SeqViennaCL;
200:   vv->ops->axpy            = VecAXPY_SeqViennaCL;
201:   vv->ops->axpby           = VecAXPBY_SeqViennaCL;
202:   vv->ops->maxpy           = VecMAXPY_SeqViennaCL;
203:   vv->ops->aypx            = VecAYPX_SeqViennaCL;
204:   vv->ops->axpbypcz        = VecAXPBYPCZ_SeqViennaCL;
205:   vv->ops->pointwisemult   = VecPointwiseMult_SeqViennaCL;
206:   vv->ops->setrandom       = VecSetRandom_SeqViennaCL;
207:   vv->ops->replacearray    = VecReplaceArray_SeqViennaCL;
208:   vv->ops->dot_local       = VecDot_SeqViennaCL;
209:   vv->ops->tdot_local      = VecTDot_SeqViennaCL;
210:   vv->ops->norm_local      = VecNorm_SeqViennaCL;
211:   vv->ops->mdot_local      = VecMDot_SeqViennaCL;
212:   vv->ops->destroy         = VecDestroy_MPIViennaCL;
213:   vv->ops->pointwisedivide = VecPointwiseDivide_SeqViennaCL;
214:   /* place array?
215:      reset array?
216:      get values?
217:   */
218:   VecViennaCLAllocateCheck(vv);
219:   vv->valid_GPU_array      = PETSC_VIENNACL_GPU;
220:   VecSet(vv,0.0);
221:   return(0);
222: }


227: PETSC_EXTERN PetscErrorCode VecCreate_ViennaCL(Vec v)
228: {
230:   PetscMPIInt    size;

233:   MPI_Comm_size(PetscObjectComm((PetscObject)v),&size);
234:   if (size == 1) {
235:     VecSetType(v,VECSEQVIENNACL);
236:   } else {
237:     VecSetType(v,VECMPIVIENNACL);
238:   }
239:   return(0);
240: }