Actual source code: textbelt.c

petsc-3.14.6 2021-03-30
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:
 23:     TextBelt is run for testing purposes only, please do not use this feature often

 25:    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
 26:    registering and using PetscTellMyCell() instead. Or email us with other alternatives we might add or make a pull request.

 28:    Developer Notes:
 29:     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
 30:        waiting for part of the message to arrive that does not exist, hence the success flg may be improperly set to false even
 31:        though the message was delivered.

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

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

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