Actual source code: const.c

petsc-3.6.1 2015-08-06
Report Typos and Errors
  2: #include <../src/vec/pf/pfimpl.h>            /*I "petscpf.h" I*/

  6: PetscErrorCode PFApply_Constant(void *value,PetscInt n,const PetscScalar *x,PetscScalar *y)
  7: {
  8:   PetscInt    i;
  9:   PetscScalar v = ((PetscScalar*)value)[0];

 12:   n *= (PetscInt) PetscRealPart(((PetscScalar*)value)[1]);
 13:   for (i=0; i<n; i++) y[i] = v;
 14:   return(0);
 15: }

 19: PetscErrorCode PFApplyVec_Constant(void *value,Vec x,Vec y)
 20: {

 24:   VecSet(y,*((PetscScalar*)value));
 25:   return(0);
 26: }
 29: PetscErrorCode PFView_Constant(void *value,PetscViewer viewer)
 30: {
 32:   PetscBool      iascii;

 35:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
 36:   if (iascii) {
 37: #if !defined(PETSC_USE_COMPLEX)
 38:     PetscViewerASCIIPrintf(viewer,"Constant = %g\n",*(double*)value);
 39: #else
 40:     PetscViewerASCIIPrintf(viewer,"Constant = %g + %gi\n",PetscRealPart(*(PetscScalar*)value),PetscImaginaryPart(*(PetscScalar*)value));
 41: #endif
 42:   }
 43:   return(0);
 44: }
 47: PetscErrorCode PFDestroy_Constant(void *value)
 48: {

 52:   PetscFree(value);
 53:   return(0);
 54: }

 58: PetscErrorCode PFSetFromOptions_Constant(PetscOptions *PetscOptionsObject,PF pf)
 59: {
 61:   PetscScalar    *value = (PetscScalar*)pf->data;

 64:   PetscOptionsHead(PetscOptionsObject,"Constant function options");
 65:   PetscOptionsScalar("-pf_constant","The constant value","None",*value,value,0);
 66:   PetscOptionsTail();
 67:   return(0);
 68: }

 72: PETSC_EXTERN PetscErrorCode PFCreate_Constant(PF pf,void *value)
 73: {
 75:   PetscScalar    *loc;

 78:   PetscMalloc1(2,&loc);
 79:   if (value) loc[0] = *(PetscScalar*)value;
 80:   else loc[0] = 0.0;
 81:   loc[1] = pf->dimout;
 82:   PFSet(pf,PFApply_Constant,PFApplyVec_Constant,PFView_Constant,PFDestroy_Constant,loc);

 84:   pf->ops->setfromoptions = PFSetFromOptions_Constant;
 85:   return(0);
 86: }

 88: /*typedef PetscErrorCode (*FCN)(void*,PetscInt,const PetscScalar*,PetscScalar*);  force argument to next function to not be extern C*/

 92: PETSC_EXTERN PetscErrorCode PFCreate_Quick(PF pf,PetscErrorCode (*function)(void*,PetscInt,const PetscScalar*,PetscScalar*))
 93: {

 97:   PFSet(pf,function,0,0,0,0);
 98:   return(0);
 99: }

101: /* -------------------------------------------------------------------------------------------------------------------*/
104: PetscErrorCode PFApply_Identity(void *value,PetscInt n,const PetscScalar *x,PetscScalar *y)
105: {
106:   PetscInt i;

109:   n *= *(PetscInt*)value;
110:   for (i=0; i<n; i++) y[i] = x[i];
111:   return(0);
112: }

116: PetscErrorCode PFApplyVec_Identity(void *value,Vec x,Vec y)
117: {

121:   VecCopy(x,y);
122:   return(0);
123: }
126: PetscErrorCode PFView_Identity(void *value,PetscViewer viewer)
127: {
129:   PetscBool      iascii;

132:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
133:   if (iascii) {
134:     PetscViewerASCIIPrintf(viewer,"Identity function\n");
135:   }
136:   return(0);
137: }
140: PetscErrorCode PFDestroy_Identity(void *value)
141: {

145:   PetscFree(value);
146:   return(0);
147: }

151: PETSC_EXTERN PetscErrorCode PFCreate_Identity(PF pf,void *value)
152: {
154:   PetscInt       *loc;

157:   if (pf->dimout != pf->dimin) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Input dimension must match output dimension for Identity function, dimin = %D dimout = %D\n",pf->dimin,pf->dimout);
158:   PetscMalloc(sizeof(PetscInt),&loc);
159:   loc[0] = pf->dimout;
160:   PFSet(pf,PFApply_Identity,PFApplyVec_Identity,PFView_Identity,PFDestroy_Identity,loc);
161:   return(0);
162: }