Actual source code: ex25.c

petsc-3.13.6 2020-09-29
Report Typos and Errors
  1: static char help[] = "Tests wrapping of math.h functions for real, complex, and scalar types \n";
  2:  #include <petscsys.h>

  4: int main(int argc,char **argv)
  5: {

  8:   PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
  9:   PetscPrintf(PETSC_COMM_WORLD,"Real tests:\n");
 10:   {
 11:     PetscReal a,b,c;
 12:     a = PetscRealConstant(0.5);
 13:     c = PetscRealConstant(2.0);

 15:     b = PetscSqrtReal(a);
 16:     PetscPrintf(PETSC_COMM_WORLD,"sqrt(%f) = %f\n",(double)a,(double)b);
 17:     b = PetscCbrtReal(a);
 18:     PetscPrintf(PETSC_COMM_WORLD,"cbrt(%f) = %f\n",(double)a,(double)b);

 20:     b = PetscHypotReal(a,c);
 21:     PetscPrintf(PETSC_COMM_WORLD,"hypot(%f,%f) = %f\n",(double)a,(double)c,(double)b);
 22:     b = PetscAtan2Real(a,c);
 23:     PetscPrintf(PETSC_COMM_WORLD,"atan2(%f,%f) = %f\n",(double)a,(double)c,(double)b);

 25:     b = PetscPowReal(a,c);
 26:     PetscPrintf(PETSC_COMM_WORLD,"pow(%f,%f) = %f\n",(double)a,(double)c,(double)b);
 27:     b = PetscExpReal(a);
 28:     PetscPrintf(PETSC_COMM_WORLD,"exp(%f) = %f\n",(double)a,(double)b);
 29:     b = PetscLogReal(a);
 30:     PetscPrintf(PETSC_COMM_WORLD,"log(%f) = %f\n",(double)a,(double)b);
 31:     b = PetscLog10Real(a);
 32:     PetscPrintf(PETSC_COMM_WORLD,"log10(%f) = %f\n",(double)a,(double)b);
 33:     b = PetscLog2Real(a);
 34:     PetscPrintf(PETSC_COMM_WORLD,"log2(%f) = %f\n",(double)a,(double)b);

 36:     b = PetscSinReal(a);
 37:     PetscPrintf(PETSC_COMM_WORLD,"sin(%f) = %f\n",(double)a,(double)b);
 38:     b = PetscCosReal(a);
 39:     PetscPrintf(PETSC_COMM_WORLD,"cos(%f) = %f\n",(double)a,(double)b);
 40:     b = PetscTanReal(a);
 41:     PetscPrintf(PETSC_COMM_WORLD,"tan(%f) = %f\n",(double)a,(double)b);

 43:     b = PetscAsinReal(a);
 44:     PetscPrintf(PETSC_COMM_WORLD,"asin(%f) = %f\n",(double)a,(double)b);
 45:     b = PetscAcosReal(a);
 46:     PetscPrintf(PETSC_COMM_WORLD,"acos(%f) = %f\n",(double)a,(double)b);
 47:     b = PetscAtanReal(a);
 48:     PetscPrintf(PETSC_COMM_WORLD,"atan(%f) = %f\n",(double)a,(double)b);

 50:     b = PetscSinhReal(a);
 51:     PetscPrintf(PETSC_COMM_WORLD,"sinh(%f) = %f\n",(double)a,(double)b);
 52:     b = PetscCoshReal(a);
 53:     PetscPrintf(PETSC_COMM_WORLD,"cosh(%f) = %f\n",(double)a,(double)b);
 54:     b = PetscTanhReal(a);
 55:     PetscPrintf(PETSC_COMM_WORLD,"tanh(%f) = %f\n",(double)a,(double)b);

 57:     b = PetscAsinhReal(a);
 58:     PetscPrintf(PETSC_COMM_WORLD,"asinh(%f) = %f\n",(double)a,(double)b);
 59:     b = PetscAcoshReal(c);
 60:     PetscPrintf(PETSC_COMM_WORLD,"acosh(%f) = %f\n",(double)c,(double)b);
 61:     b = PetscAtanhReal(a);
 62:     PetscPrintf(PETSC_COMM_WORLD,"atanh(%f) = %f\n",(double)a,(double)b);

 64:     b = PetscCeilReal(a);
 65:     PetscPrintf(PETSC_COMM_WORLD,"ceil(%f) = %f\n",(double)a,(double)b);
 66:     b = PetscFloorReal(a);
 67:     PetscPrintf(PETSC_COMM_WORLD,"floor(%f) = %f\n",(double)a,(double)b);
 68:     b = PetscFmodReal(a,c);
 69:     PetscPrintf(PETSC_COMM_WORLD,"fmod(%f,%f) = %f\n",(double)a,(double)c,(double)b);
 70:   }
 71:   PetscPrintf(PETSC_COMM_WORLD,"Scalar tests:\n");
 72:   {
 73:     PetscScalar a,b,c;
 74:     a = PetscRealConstant(0.5);
 75:     c = PetscRealConstant(2.0);

 77:     b = PetscAbsScalar(a);
 78:     PetscPrintf(PETSC_COMM_WORLD,"abs(%f) = %f\n",(double)PetscRealPart(a),(double)PetscRealPart(b));
 79:     b = PetscArgScalar(a);
 80:     PetscPrintf(PETSC_COMM_WORLD,"arg(%f) = %f\n",(double)PetscRealPart(a),(double)PetscRealPart(b));
 81:     b = PetscConj(a);
 82:     PetscPrintf(PETSC_COMM_WORLD,"conj(%f) = %f\n",(double)PetscRealPart(a),(double)PetscRealPart(b));

 84:     b = PetscSqrtScalar(a);
 85:     PetscPrintf(PETSC_COMM_WORLD,"sqrt(%f) = %f\n",(double)PetscRealPart(a),(double)PetscRealPart(b));

 87:     b = PetscPowScalar(a,c);
 88:     PetscPrintf(PETSC_COMM_WORLD,"pow(%f,%f) = %f\n",(double)PetscRealPart(a),(double)PetscRealPart(c),(double)PetscRealPart(b));
 89:     b = PetscExpScalar(a);
 90:     PetscPrintf(PETSC_COMM_WORLD,"exp(%f) = %f\n",(double)PetscRealPart(a),(double)PetscRealPart(b));
 91:     b = PetscLogScalar(a);
 92:     PetscPrintf(PETSC_COMM_WORLD,"log(%f) = %f\n",(double)PetscRealPart(a),(double)PetscRealPart(b));

 94:     b = PetscSinScalar(a);
 95:     PetscPrintf(PETSC_COMM_WORLD,"sin(%f) = %f\n",(double)PetscRealPart(a),(double)PetscRealPart(b));
 96:     b = PetscCosScalar(a);
 97:     PetscPrintf(PETSC_COMM_WORLD,"cos(%f) = %f\n",(double)PetscRealPart(a),(double)PetscRealPart(b));
 98:     b = PetscTanScalar(a);
 99:     PetscPrintf(PETSC_COMM_WORLD,"tan(%f) = %f\n",(double)PetscRealPart(a),(double)PetscRealPart(b));

101:     b = PetscAsinScalar(a);
102:     PetscPrintf(PETSC_COMM_WORLD,"asin(%f) = %f\n",(double)PetscRealPart(a),(double)PetscRealPart(b));
103:     b = PetscAcosScalar(a);
104:     PetscPrintf(PETSC_COMM_WORLD,"acos(%f) = %f\n",(double)PetscRealPart(a),(double)PetscRealPart(b));
105:     b = PetscAtanScalar(a);
106:     PetscPrintf(PETSC_COMM_WORLD,"atan(%f) = %f\n",(double)PetscRealPart(a),(double)PetscRealPart(b));

108:     b = PetscSinhScalar(a);
109:     PetscPrintf(PETSC_COMM_WORLD,"sinh(%f) = %f\n",(double)PetscRealPart(a),(double)PetscRealPart(b));
110:     b = PetscCoshScalar(a);
111:     PetscPrintf(PETSC_COMM_WORLD,"cosh(%f) = %f\n",(double)PetscRealPart(a),(double)PetscRealPart(b));
112:     b = PetscTanhScalar(a);
113:     PetscPrintf(PETSC_COMM_WORLD,"tanh(%f) = %f\n",(double)PetscRealPart(a),(double)PetscRealPart(b));

115:     b = PetscAsinhScalar(a);
116:     PetscPrintf(PETSC_COMM_WORLD,"asinh(%f) = %f\n",(double)PetscRealPart(a),(double)PetscRealPart(b));
117:     b = PetscAcoshScalar(c);
118:     PetscPrintf(PETSC_COMM_WORLD,"acosh(%f) = %f\n",(double)PetscRealPart(c),(double)PetscRealPart(b));
119:     b = PetscAtanhScalar(a);
120:     PetscPrintf(PETSC_COMM_WORLD,"atanh(%f) = %f\n",(double)PetscRealPart(a),(double)PetscRealPart(b));
121:   }
122:   PetscFinalize();
123:   return ierr;
124: }


127: /*TEST

129:    test:

131: TEST*/