Actual source code: textbelt.c

petsc-3.8.4 2018-03-24
Report Typos and Errors

  2:  #include <petscwebclient.h>

  4: /*@C
  5:      PetscTextBelt - Sends an SMS to an American/Canadian phone number

  7:    Not collective, only the first process in MPI_Comm does anything

  9:    Input Parameters:
 10: +  comm - the MPI communicator
 11: .  number - the 10 digit telephone number
 12: -  message - the message

 14:    Output Parameter:
 15: .   flg - PETSC_TRUE if the text was sent

 17:    Options Database:
 18: .   -textbelt <phonenumber[,message]>

 20:    Level: intermediate

 22:    Notes: TextBelt is run for testing purposes only, please do not use this feature often

 24:    As of November 2016 this service does not seem to be actually transmitting the SMS, which is unfortunate since it is such a great service. Consider
 25:    registering and using PetscTellMyCell() instead. Or email us with other alternatives we might add or make a pull request.

 27:    Developer Notes:  I do not know how to make the buff[] long enough to receive the "success" string but short enough that the code does not hang
 28:        waiting for part of the message to arrive that does not exist, hence the success flg may be improperly set to false even
 29:        though the message was delivered.

 31: .seealso: PetscTellMyCell(), PetscOpenSocket(), PetscHTTPRequest()
 32: @*/
 33: PetscErrorCode PetscTextBelt(MPI_Comm comm,const char number[],const char message[],PetscBool *flg)
 34: {
 36:   size_t         nlen,mlen,blen;
 37:   PetscMPIInt    rank;

 40:   PetscStrlen(number,&nlen);
 41:   if (nlen != 10) SETERRQ1(comm,PETSC_ERR_ARG_WRONG,"Number %s is not ten digits",number);
 42:   PetscStrlen(message,&mlen);
 43:   if (mlen > 100) SETERRQ1(comm,PETSC_ERR_ARG_WRONG,"Message  %s is too long",message);
 44:   MPI_Comm_rank(comm,&rank);
 45:   if (!rank) {
 46:     int       sock;
 47:     char      buff[474],*body;
 48:     PetscInt  i;

 50:     PetscMalloc1(mlen+nlen+100,&body);
 51:     PetscStrcpy(body,"number=");
 52:     PetscStrcat(body,number);
 53:     PetscStrcat(body,"&");
 54:     PetscStrcat(body,"message=");
 55:     PetscStrcat(body,message);
 56:     PetscStrlen(body,&blen);
 57:     for (i=0; i<(int)blen; i++) {
 58:       if (body[i] == ' ') body[i] = '+';
 59:     }
 60:     PetscOpenSocket("textbelt.com",80,&sock);
 61:     PetscHTTPRequest("POST","textbelt.com/text",NULL,"application/x-www-form-urlencoded",body,sock,buff,sizeof(buff));
 62:     close(sock);
 63:     PetscFree(body);
 64:     if (flg) {
 65:       char *found;
 66:       PetscStrstr(buff,"\"success\":tr",&found);
 67:       *flg = found ? PETSC_TRUE : PETSC_FALSE;
 68:     }
 69:   }
 70:   return(0);
 71: }