Actual source code: textbelt.c
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 Key:
18: . -textbelt <phonenumber[,message]> - sends a message to this number when the program ends
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 Note:
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: {
37: size_t nlen, mlen, blen;
38: PetscMPIInt rank;
40: PetscFunctionBegin;
41: PetscCall(PetscStrlen(number, &nlen));
42: PetscCheck(nlen == 10, comm, PETSC_ERR_ARG_WRONG, "Number %s is not ten digits", number);
43: PetscCall(PetscStrlen(message, &mlen));
44: PetscCheck(mlen <= 100, comm, PETSC_ERR_ARG_WRONG, "Message %s is too long", message);
45: PetscCallMPI(MPI_Comm_rank(comm, &rank));
46: if (rank == 0) {
47: int sock;
48: char buff[474], *body;
49: PetscInt i;
51: blen = mlen + nlen + 100;
52: PetscCall(PetscMalloc1(blen, &body));
53: PetscCall(PetscStrncpy(body, "number=", blen));
54: PetscCall(PetscStrlcat(body, number, blen));
55: PetscCall(PetscStrlcat(body, "&", blen));
56: PetscCall(PetscStrlcat(body, "message=", blen));
57: PetscCall(PetscStrlcat(body, message, blen));
58: PetscCall(PetscStrlen(body, &blen));
59: for (i = 0; i < (int)blen; i++) {
60: if (body[i] == ' ') body[i] = '+';
61: }
62: PetscCall(PetscOpenSocket("textbelt.com", 80, &sock));
63: PetscCall(PetscHTTPRequest("POST", "textbelt.com/text", NULL, "application/x-www-form-urlencoded", body, sock, buff, sizeof(buff)));
64: close(sock);
65: PetscCall(PetscFree(body));
66: if (flg) {
67: char *found;
68: PetscCall(PetscStrstr(buff, "\"success\":tr", &found));
69: *flg = found ? PETSC_TRUE : PETSC_FALSE;
70: }
71: }
72: PetscFunctionReturn(PETSC_SUCCESS);
73: }