Actual source code: textbelt.c

petsc-3.7.3 2016-08-01
Report Typos and Errors
  2: #include <petscwebclient.h>

  6: /*@C
  7:      PetscTextBelt - Sends an SMS to an American/Canadian phone number

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

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

 16:    Output Parameter:
 17: .   flg - PETSC_TRUE if the text was sent

 19:    Level: intermediate

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

 23:    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
 24:        waiting for part of the message to arrive that does not exist, hence the success flg may be improperly set to false even
 25:        though the message was delivered.

 27: .seealso: PetscOpenSocket(), PetscHTTPRequest()
 28: @*/
 29: PetscErrorCode PetscTextBelt(MPI_Comm comm,const char number[],const char message[],PetscBool *flg)
 30: {
 32:   size_t         nlen,mlen,blen;
 33:   PetscMPIInt    rank;

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

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