Actual source code: const.c

petsc-3.11.4 2019-09-28
Report Typos and Errors

  2:  #include <../src/vec/pf/pfimpl.h>

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

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

 15: static PetscErrorCode PFApplyVec_Constant(void *value,Vec x,Vec y)
 16: {

 20:   VecSet(y,*((PetscScalar*)value));
 21:   return(0);
 22: }
 23: PetscErrorCode PFView_Constant(void *value,PetscViewer viewer)
 24: {
 26:   PetscBool      iascii;

 29:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
 30:   if (iascii) {
 31: #if !defined(PETSC_USE_COMPLEX)
 32:     PetscViewerASCIIPrintf(viewer,"Constant = %g\n",*(double*)value);
 33: #else
 34:     PetscViewerASCIIPrintf(viewer,"Constant = %g + %gi\n",PetscRealPart(*(PetscScalar*)value),PetscImaginaryPart(*(PetscScalar*)value));
 35: #endif
 36:   }
 37:   return(0);
 38: }
 39: static PetscErrorCode PFDestroy_Constant(void *value)
 40: {

 44:   PetscFree(value);
 45:   return(0);
 46: }

 48: static PetscErrorCode PFSetFromOptions_Constant(PetscOptionItems *PetscOptionsObject,PF pf)
 49: {
 51:   PetscScalar    *value = (PetscScalar*)pf->data;

 54:   PetscOptionsHead(PetscOptionsObject,"Constant function options");
 55:   PetscOptionsScalar("-pf_constant","The constant value","None",*value,value,0);
 56:   PetscOptionsTail();
 57:   return(0);
 58: }

 60: PETSC_EXTERN PetscErrorCode PFCreate_Constant(PF pf,void *value)
 61: {
 63:   PetscScalar    *loc;

 66:   PetscMalloc1(2,&loc);
 67:   if (value) loc[0] = *(PetscScalar*)value;
 68:   else loc[0] = 0.0;
 69:   loc[1] = pf->dimout;
 70:   PFSet(pf,PFApply_Constant,PFApplyVec_Constant,PFView_Constant,PFDestroy_Constant,loc);

 72:   pf->ops->setfromoptions = PFSetFromOptions_Constant;
 73:   return(0);
 74: }

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

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

 83:   PFSet(pf,function,0,0,0,0);
 84:   return(0);
 85: }

 87: /* -------------------------------------------------------------------------------------------------------------------*/
 88: static PetscErrorCode PFApply_Identity(void *value,PetscInt n,const PetscScalar *x,PetscScalar *y)
 89: {
 90:   PetscInt i;

 93:   n *= *(PetscInt*)value;
 94:   for (i=0; i<n; i++) y[i] = x[i];
 95:   return(0);
 96: }

 98: static PetscErrorCode PFApplyVec_Identity(void *value,Vec x,Vec y)
 99: {

103:   VecCopy(x,y);
104:   return(0);
105: }
106: static PetscErrorCode PFView_Identity(void *value,PetscViewer viewer)
107: {
109:   PetscBool      iascii;

112:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
113:   if (iascii) {
114:     PetscViewerASCIIPrintf(viewer,"Identity function\n");
115:   }
116:   return(0);
117: }
118: static PetscErrorCode PFDestroy_Identity(void *value)
119: {

123:   PetscFree(value);
124:   return(0);
125: }

127: PETSC_EXTERN PetscErrorCode PFCreate_Identity(PF pf,void *value)
128: {
130:   PetscInt       *loc;

133:   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);
134:   PetscNew(&loc);
135:   loc[0] = pf->dimout;
136:   PFSet(pf,PFApply_Identity,PFApplyVec_Identity,PFView_Identity,PFDestroy_Identity,loc);
137:   return(0);
138: }